]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PHOS/AliPHOSClusterizerv1.cxx
Updated Linkdef and libTOF.pkg
[u/mrichter/AliRoot.git] / PHOS / AliPHOSClusterizerv1.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /* $Log:
19    1 October 2000. Yuri Kharlov:
20      AreNeighbours()
21      PPSD upper layer is considered if number of layers>1
22
23    18 October 2000. Yuri Kharlov:
24      AliPHOSClusterizerv1()
25      CPV clusterizing parameters added
26
27      MakeClusters()
28      After first PPSD digit remove EMC digits only once
29 */
30 //*-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
31 //////////////////////////////////////////////////////////////////////////////
32 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
33 //  unfolds the clusters having several local maxima.  
34 //  Results are stored in TreeR#, branches PHOSEmcRP (EMC recPoints),
35 //  PHOSCpvRP (CPV RecPoints) and AliPHOSClusterizer (Clusterizer with all 
36 //  parameters including input digits branch title, thresholds etc.)
37 //  This TTask is normally called from Reconstructioner, but can as well be used in 
38 //  standalone mode.
39 // Use Case:
40 //  root [0] AliPHOSClusterizerv1 * cl = new AliPHOSClusterizerv1("galice.root")  
41 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
42 //               //reads gAlice from header file "..."                      
43 //  root [1] cl->ExecuteTask()  
44 //               //finds RecPoints in all events stored in galice.root
45 //  root [2] cl->SetDigitsBranch("digits2") 
46 //               //sets another title for Digitis (input) branch
47 //  root [3] cl->SetRecPointsBranch("recp2")  
48 //               //sets another title four output branches
49 //  root [4] cl->SetEmcLocalMaxCut(0.03)  
50 //               //set clusterization parameters
51 //  root [5] cl->ExecuteTask("deb all time")  
52 //               //once more finds RecPoints options are 
53 //               // deb - print number of found rec points
54 //               // deb all - print number of found RecPoints and some their characteristics 
55 //               // time - print benchmarking results
56
57 // --- ROOT system ---
58
59 #include "TROOT.h" 
60 #include "TFile.h" 
61 #include "TFolder.h" 
62 #include "TMath.h" 
63 #include "TMinuit.h"
64 #include "TTree.h" 
65 #include "TSystem.h" 
66 #include "TBenchmark.h"
67
68 // --- Standard library ---
69
70 #include <iostream.h>
71 #include <iomanip.h>
72
73 // --- AliRoot header files ---
74
75 #include "AliPHOSClusterizerv1.h"
76 #include "AliPHOSCpvRecPoint.h"
77 #include "AliPHOSDigit.h"
78 #include "AliPHOSDigitizer.h"
79 #include "AliPHOSEmcRecPoint.h"
80 #include "AliPHOS.h"
81 #include "AliPHOSGetter.h"
82 #include "AliRun.h"
83
84 ClassImp(AliPHOSClusterizerv1)
85   
86 //____________________________________________________________________________
87   AliPHOSClusterizerv1::AliPHOSClusterizerv1() : AliPHOSClusterizer()
88 {
89   // default ctor (to be used mainly by Streamer)
90   
91   fNumberOfCpvClusters     = 0 ; 
92   fNumberOfEmcClusters     = 0 ; 
93   
94   fCpvClusteringThreshold  = 0.0;
95   fEmcClusteringThreshold  = 0.2;   
96   
97   fEmcLocMaxCut            = 0.03 ;
98   fCpvLocMaxCut            = 0.03 ;
99   
100   fW0                      = 4.5 ;
101   fW0CPV                   = 4.0 ;
102
103   fEmcTimeGate             = 1.e-8 ; 
104
105   fToUnfold                = kTRUE ;
106
107   fHeaderFileName          = "" ; 
108   fDigitsBranchTitle       = "" ;
109   fRecPointsInRun          = 0 ;   
110 }
111
112 //____________________________________________________________________________
113 AliPHOSClusterizerv1::AliPHOSClusterizerv1(const char* headerFile,const char* name)
114 :AliPHOSClusterizer(headerFile, name)
115 {
116   // ctor with the indication of the file where header Tree and digits Tree are stored
117   
118
119   fNumberOfCpvClusters     = 0 ; 
120   fNumberOfEmcClusters     = 0 ; 
121   
122   fCpvClusteringThreshold  = 0.0;
123   fEmcClusteringThreshold  = 0.2;   
124   
125   fEmcLocMaxCut            = 0.03 ;
126   fCpvLocMaxCut            = 0.03 ;
127   
128   fW0                      = 4.5 ;
129   fW0CPV                   = 4.0 ;
130
131   fEmcTimeGate             = 1.e-8 ; 
132   
133   fToUnfold                = kTRUE ;
134   
135   fHeaderFileName          = GetTitle() ; 
136   fDigitsBranchTitle       = GetName() ;
137   
138   TString clusterizerName( GetName()) ; 
139   clusterizerName.Append(":") ; 
140   clusterizerName.Append(Version()) ; 
141   SetName(clusterizerName) ;
142   fRecPointsInRun          = 0 ; 
143
144   Init() ;
145
146 }
147 //____________________________________________________________________________
148   AliPHOSClusterizerv1::~AliPHOSClusterizerv1()
149 {
150 }
151 //____________________________________________________________________________
152 Float_t  AliPHOSClusterizerv1::Calibrate(Int_t amp, Int_t absId) const
153 {
154   if(absId <= fEmcCrystals) //calibrate as EMC 
155     return fADCpedestalEmc + amp*fADCchanelEmc ;       
156   else //Digitize as CPV
157     return fADCpedestalCpv+ amp*fADCchanelCpv ;       
158 }
159 //____________________________________________________________________________
160 void AliPHOSClusterizerv1::Exec(Option_t * option)
161 {
162   // Steering method
163
164   if( strcmp(GetName(), "")== 0 ) 
165     Init() ;
166
167   if(strstr(option,"tim"))
168     gBenchmark->Start("PHOSClusterizer"); 
169   
170   if(strstr(option,"print"))
171     Print("") ; 
172
173   gAlice->GetEvent(0) ;
174   
175  //check, if the branch with name of this" already exits?
176   TObjArray * lob = (TObjArray*)gAlice->TreeR()->GetListOfBranches() ;
177   TIter next(lob) ; 
178   TBranch * branch = 0 ;  
179   Bool_t phosemcfound = kFALSE, phoscpvfound = kFALSE, clusterizerfound = kFALSE ; 
180
181   TString branchname = GetName() ;
182   branchname.Remove(branchname.Index(Version())-1) ;
183   
184   while ( (branch = (TBranch*)next()) && (!phosemcfound || !phoscpvfound || !clusterizerfound) ) {
185     if ( (strcmp(branch->GetName(), "PHOSEmcRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
186       phosemcfound = kTRUE ;
187     
188     else if ( (strcmp(branch->GetName(), "PHOSCpvRP")==0) && (strcmp(branch->GetTitle(), branchname.Data())==0) ) 
189       phoscpvfound = kTRUE ;
190    
191     else if ((strcmp(branch->GetName(), "AliPHOSClusterizer")==0) && (strcmp(branch->GetTitle(), GetName())==0) ) 
192       clusterizerfound = kTRUE ; 
193   }
194
195   if ( phoscpvfound || phosemcfound || clusterizerfound ) {
196     cerr << "WARNING: AliPHOSClusterizer::Exec -> Emc(Cpv)RecPoints and/or Clusterizer branch with name " 
197          << branchname.Data() << " already exits" << endl ;
198     return ; 
199   }       
200
201   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
202   Int_t nevents = (Int_t) gAlice->TreeE()->GetEntries() ;
203   Int_t ievent ;
204
205   for(ievent = 0; ievent < nevents; ievent++){
206
207     if(ievent == 0) GetCalibrationParameters() ;
208
209     fNumberOfEmcClusters  = 0 ;
210     fNumberOfCpvClusters  = 0 ;
211    
212     gime->Event(ievent,"D") ;
213     
214     //    if(!ReadDigits(ievent))   continue;  //reads digits for event ievent
215     
216     MakeClusters() ;
217     
218     if(fToUnfold)             MakeUnfolding() ;
219
220     WriteRecPoints(ievent) ;
221
222     if(strstr(option,"deb"))  
223       PrintRecPoints(option) ;
224
225     //increment the total number of digits per run 
226     fRecPointsInRun += gime->EmcRecPoints()->GetEntriesFast() ;  
227     fRecPointsInRun += gime->CpvRecPoints()->GetEntriesFast() ;  
228  }
229   
230   if(strstr(option,"tim")){
231     gBenchmark->Stop("PHOSClusterizer");
232     cout << "AliPHOSClusterizer:" << endl ;
233     cout << "  took " << gBenchmark->GetCpuTime("PHOSClusterizer") << " seconds for Clusterizing " 
234          <<  gBenchmark->GetCpuTime("PHOSClusterizer")/nevents << " seconds per event " << endl ;
235     cout << endl ;
236   }
237   
238 }
239
240 //____________________________________________________________________________
241 Bool_t AliPHOSClusterizerv1::FindFit(AliPHOSEmcRecPoint * emcRP, int * maxAt, Float_t * maxAtEnergy,
242                                     Int_t nPar, Float_t * fitparameters) const
243
244   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
245   // The initial values for fitting procedure are set equal to the positions of local maxima.
246   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
247
248   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
249   TClonesArray * digits = gime->Digits() ; 
250   
251
252   gMinuit->mncler();                     // Reset Minuit's list of paramters
253   gMinuit->SetPrintLevel(-1) ;           // No Printout
254   gMinuit->SetFCN(AliPHOSClusterizerv1::UnfoldingChiSquare) ;  
255                                          // To set the address of the minimization function 
256
257   TList * toMinuit = new TList();
258   toMinuit->AddAt(emcRP,0) ;
259   toMinuit->AddAt(digits,1) ;
260   
261   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
262
263   // filling initial values for fit parameters
264   AliPHOSDigit * digit ;
265
266   Int_t ierflg  = 0; 
267   Int_t index   = 0 ;
268   Int_t nDigits = (Int_t) nPar / 3 ;
269
270   Int_t iDigit ;
271
272   const AliPHOSGeometry * geom = gime->PHOSGeometry() ; 
273
274   for(iDigit = 0; iDigit < nDigits; iDigit++){
275     digit = (AliPHOSDigit *) maxAt[iDigit]; 
276
277     Int_t relid[4] ;
278     Float_t x = 0.;
279     Float_t z = 0.;
280     geom->AbsToRelNumbering(digit->GetId(), relid) ;
281     geom->RelPosInModule(relid, x, z) ;
282
283     Float_t energy = maxAtEnergy[iDigit] ;
284
285     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
286     index++ ;   
287     if(ierflg != 0){ 
288       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : x = " << x << endl ;
289       return kFALSE;
290     }
291     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
292     index++ ;   
293     if(ierflg != 0){
294       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : z = " << z << endl ;
295       return kFALSE;
296     }
297     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
298     index++ ;   
299     if(ierflg != 0){
300       cout << "PHOS Unfolding>  Unable to set initial value for fit procedure : energy = " << energy << endl ;      
301       return kFALSE;
302     }
303   }
304
305   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
306                       //  depends on it. 
307   Double_t p1 = 1.0 ;
308   Double_t p2 = 0.0 ;
309
310   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
311   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
312   gMinuit->SetMaxIterations(5);
313   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
314
315   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
316
317   if(ierflg == 4){  // Minimum not found   
318     cout << "PHOS Unfolding>  Fit not converged, cluster abandoned "<< endl ;      
319     return kFALSE ;
320   }            
321   for(index = 0; index < nPar; index++){
322     Double_t err ;
323     Double_t val ;
324     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
325     fitparameters[index] = val ;
326    }
327
328   delete toMinuit ;
329   return kTRUE;
330
331 }
332
333 //____________________________________________________________________________
334 void AliPHOSClusterizerv1::GetCalibrationParameters() 
335 {
336   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ;
337   TString branchname = GetName() ;
338   branchname.Remove(branchname.Index(Version())-1) ;
339   AliPHOSDigitizer * dig = gime->Digitizer(branchname) ;
340   fADCchanelEmc   = dig->GetEMCchannel() ;
341   fADCpedestalEmc = dig->GetEMCpedestal();
342
343   fADCchanelCpv   = dig->GetCPVchannel() ;
344   fADCpedestalCpv = dig->GetCPVpedestal() ; 
345
346 }
347 //____________________________________________________________________________
348 void AliPHOSClusterizerv1::Init()
349 {
350   // Make all memory allocations which can not be done in default constructor.
351   // Attach the Clusterizer task to the list of PHOS tasks
352   
353   if ( strcmp(GetTitle(), "") == 0 )
354     SetTitle("galice.root") ;
355
356   TString branchname = GetName() ;
357   branchname.Remove(branchname.Index(Version())-1) ;
358
359   AliPHOSGetter * gime = AliPHOSGetter::GetInstance(GetTitle(), branchname) ; 
360   if ( gime == 0 ) {
361     cerr << "ERROR: AliPHOSClusterizerv1::Init -> Could not obtain the Getter object !" << endl ; 
362     return ;
363   } 
364     
365   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
366   fEmcCrystals = geom->GetNModules() *  geom->GetNCristalsInModule() ;
367
368   if(!gMinuit) 
369     gMinuit = new TMinuit(100) ;
370
371   gime->PostClusterizer(this) ;
372   // create a folder on the white board //YSAlice/WhiteBoard/RecPoints/PHOS/recpointsName
373   gime->PostRecPoints(branchname ) ;
374
375   gime->PostDigits(branchname) ;
376   gime->PostDigitizer(branchname) ;
377   
378 }
379
380 //____________________________________________________________________________
381 Int_t AliPHOSClusterizerv1::AreNeighbours(AliPHOSDigit * d1, AliPHOSDigit * d2)const
382 {
383   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
384   //                                       = 1 are neighbour
385   //                                       = 2 are not neighbour but do not continue searching
386   // neighbours are defined as digits having at least a common vertex 
387   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
388   //                                      which is compared to a digit (d2)  not yet in a cluster  
389
390   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
391
392   Int_t rv = 0 ; 
393
394   Int_t relid1[4] ; 
395   geom->AbsToRelNumbering(d1->GetId(), relid1) ; 
396
397   Int_t relid2[4] ; 
398   geom->AbsToRelNumbering(d2->GetId(), relid2) ; 
399  
400   if ( (relid1[0] == relid2[0]) && (relid1[1]==relid2[1]) ) { // inside the same PHOS module 
401     Int_t rowdiff = TMath::Abs( relid1[2] - relid2[2] ) ;  
402     Int_t coldiff = TMath::Abs( relid1[3] - relid2[3] ) ;  
403     
404     if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
405       if((relid1[1] != 0) || (TMath::Abs(d1->GetTime() - d2->GetTime() ) < fEmcTimeGate))
406       rv = 1 ; 
407     }
408     else {
409       if((relid2[2] > relid1[2]) && (relid2[3] > relid1[3]+1)) 
410         rv = 2; //  Difference in row numbers is too large to look further 
411     }
412
413   } 
414   else {
415     
416     if( (relid1[0] < relid2[0]) || (relid1[1] != relid2[1]) )  
417       rv=2 ;
418
419   }
420
421   return rv ; 
422 }
423
424
425 //____________________________________________________________________________
426 Bool_t AliPHOSClusterizerv1::IsInEmc(AliPHOSDigit * digit) const
427 {
428   // Tells if (true) or not (false) the digit is in a PHOS-EMC module
429  
430   Bool_t rv = kFALSE ; 
431   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
432
433   Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();  
434
435   if(digit->GetId() <= nEMC )   rv = kTRUE; 
436
437   return rv ; 
438 }
439
440 //____________________________________________________________________________
441 Bool_t AliPHOSClusterizerv1::IsInCpv(AliPHOSDigit * digit) const
442 {
443   // Tells if (true) or not (false) the digit is in a PHOS-CPV module
444  
445   Bool_t rv = kFALSE ; 
446   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ;
447
448   Int_t nEMC = geom->GetNModules()*geom->GetNPhi()*geom->GetNZ();
449
450   if(digit->GetId() > nEMC )   rv = kTRUE;
451
452   return rv ; 
453 }
454
455 //____________________________________________________________________________
456 void AliPHOSClusterizerv1::WriteRecPoints(Int_t event)
457 {
458
459   // Creates new branches with given title
460   // fills and writes into TreeR.
461
462   TString branchName(GetName() ) ;
463   branchName.Remove(branchName.Index(Version())-1) ;
464
465   AliPHOSGetter *gime = AliPHOSGetter::GetInstance() ; 
466   TObjArray * emcRecPoints = gime->EmcRecPoints(branchName) ; 
467   TObjArray * cpvRecPoints = gime->CpvRecPoints(branchName) ; 
468   TClonesArray * digits = gime->Digits(branchName) ; 
469
470   Int_t index ;
471   //Evaluate position, dispersion and other RecPoint properties...
472   for(index = 0; index < emcRecPoints->GetEntries(); index++)
473     ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->EvalAll(fW0,digits) ;
474
475   emcRecPoints->Sort() ;
476
477   for(index = 0; index < emcRecPoints->GetEntries(); index++)
478     ((AliPHOSEmcRecPoint *)emcRecPoints->At(index))->SetIndexInList(index) ;
479
480   emcRecPoints->Expand(emcRecPoints->GetEntriesFast()) ; 
481
482   //Now the same for CPV
483   for(index = 0; index < cpvRecPoints->GetEntries(); index++)
484     ((AliPHOSRecPoint *)cpvRecPoints->At(index))->EvalAll(fW0CPV,digits)  ;
485
486   cpvRecPoints->Sort() ;
487
488   for(index = 0; index < cpvRecPoints->GetEntries(); index++)
489     ((AliPHOSRecPoint *)cpvRecPoints->At(index))->SetIndexInList(index) ;
490
491   cpvRecPoints->Expand(cpvRecPoints->GetEntriesFast()) ;
492   
493   //Make branches in TreeR for RecPoints and Clusterizer
494   char * filename = 0;
495   if(gSystem->Getenv("CONFIG_SPLIT_FILE")!=0){   //generating file name
496     filename = new char[strlen(gAlice->GetBaseFile())+20] ;
497     sprintf(filename,"%s/PHOS.Reco.root",gAlice->GetBaseFile()) ;
498   }
499   
500   //Make new branches
501   TDirectory *cwd = gDirectory;
502   
503  
504   Int_t bufferSize = 32000 ;    
505   Int_t splitlevel = 0 ;
506
507   //First EMC
508   TBranch * emcBranch = gAlice->TreeR()->Branch("PHOSEmcRP","TObjArray",&emcRecPoints,bufferSize,splitlevel);
509   emcBranch->SetTitle(branchName);
510   if (filename) {
511     emcBranch->SetFile(filename);
512     TIter next( emcBranch->GetListOfBranches());
513     TBranch * sb ;
514     while ((sb=(TBranch*)next())) {
515       sb->SetFile(filename);
516     }   
517     
518     cwd->cd();
519   }
520     
521   //Now CPV branch
522   TBranch * cpvBranch = gAlice->TreeR()->Branch("PHOSCpvRP","TObjArray",&cpvRecPoints,bufferSize,splitlevel);
523   cpvBranch->SetTitle(branchName);
524   if (filename) {
525     cpvBranch->SetFile(filename);
526     TIter next( cpvBranch->GetListOfBranches());
527     TBranch * sb;
528     while ((sb=(TBranch*)next())) {
529       sb->SetFile(filename);
530     }   
531     cwd->cd();
532   } 
533     
534   //And Finally  clusterizer branch
535   AliPHOSClusterizerv1 * cl = (AliPHOSClusterizerv1*)gime->Clusterizer(branchName) ;
536   TBranch * clusterizerBranch = gAlice->TreeR()->Branch("AliPHOSClusterizer","AliPHOSClusterizerv1",
537                                               &cl,bufferSize,splitlevel);
538   clusterizerBranch->SetTitle(branchName);
539   if (filename) {
540     clusterizerBranch->SetFile(filename);
541     TIter next( clusterizerBranch->GetListOfBranches());
542     TBranch * sb ;
543     while ((sb=(TBranch*)next())) {
544       sb->SetFile(filename);
545     }   
546     cwd->cd();
547   }
548   emcBranch        ->Fill() ;
549   cpvBranch        ->Fill() ;
550   clusterizerBranch->Fill() ;
551
552   gAlice->TreeR()->Write(0,kOverwrite) ;  
553   
554 }
555
556 //____________________________________________________________________________
557 void AliPHOSClusterizerv1::MakeClusters()
558 {
559   // Steering method to construct the clusters stored in a list of Reconstructed Points
560   // A cluster is defined as a list of neighbour digits
561   TString branchName(GetName()) ; 
562   branchName.Remove(branchName.Index(Version())-1) ; 
563
564   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
565
566   TObjArray * emcRecPoints = gime->EmcRecPoints(branchName) ; 
567   TObjArray * cpvRecPoints = gime->CpvRecPoints(branchName) ; 
568   emcRecPoints->Delete() ;
569   cpvRecPoints->Delete() ;
570   
571   TClonesArray * digits = gime->Digits(branchName) ; 
572   TClonesArray * digitsC =  (TClonesArray*)digits->Clone() ;
573   
574  
575   // Clusterization starts  
576
577   TIter nextdigit(digitsC) ; 
578   AliPHOSDigit * digit ; 
579   Bool_t notremoved = kTRUE ;
580
581   while ( (digit = (AliPHOSDigit *)nextdigit()) ) { // scan over the list of digitsC
582     AliPHOSRecPoint * clu = 0 ; 
583
584     TArrayI clusterdigitslist(1500) ;   
585     Int_t index ;
586
587     if (( IsInEmc (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fEmcClusteringThreshold  ) || 
588         ( IsInCpv (digit) && Calibrate(digit->GetAmp(),digit->GetId()) > fCpvClusteringThreshold  ) ) {
589          
590       Int_t iDigitInCluster = 0 ; 
591       
592       if  ( IsInEmc(digit) ) {   
593         // start a new EMC RecPoint
594         if(fNumberOfEmcClusters >= emcRecPoints->GetSize()) 
595           emcRecPoints->Expand(2*fNumberOfEmcClusters+1) ;
596         
597         emcRecPoints->AddAt(new  AliPHOSEmcRecPoint(""), fNumberOfEmcClusters) ;
598         clu = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters) ; 
599         fNumberOfEmcClusters++ ; 
600         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId())) ; 
601         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList() ;  
602         iDigitInCluster++ ; 
603         digitsC->Remove(digit) ; 
604
605       } else { 
606         
607         // start a new CPV cluster
608         if(fNumberOfCpvClusters >= cpvRecPoints->GetSize()) 
609           cpvRecPoints->Expand(2*fNumberOfCpvClusters+1);
610
611         cpvRecPoints->AddAt(new AliPHOSCpvRecPoint(""), fNumberOfCpvClusters) ;
612
613         clu =  (AliPHOSCpvRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters)  ;  
614         fNumberOfCpvClusters++ ; 
615         clu->AddDigit(*digit, Calibrate(digit->GetAmp(),digit->GetId()) ) ;     
616         clusterdigitslist[iDigitInCluster] = digit->GetIndexInList()  ; 
617         iDigitInCluster++ ; 
618         digitsC->Remove(digit) ; 
619         nextdigit.Reset() ;
620         
621         // Here we remove remaining EMC digits, which cannot make a cluster
622         
623         if( notremoved ) { 
624           while( ( digit = (AliPHOSDigit *)nextdigit() ) ) {
625             if( IsInEmc(digit) ) 
626               digitsC->Remove(digit) ;
627             else 
628               break ;
629           }
630           notremoved = kFALSE ;
631         }
632         
633       } // else        
634       
635       nextdigit.Reset() ;
636       
637       AliPHOSDigit * digitN ; 
638       index = 0 ;
639       while (index < iDigitInCluster){ // scan over digits already in cluster 
640         digit =  (AliPHOSDigit*)digits->At(clusterdigitslist[index])  ;      
641         index++ ; 
642         while ( (digitN = (AliPHOSDigit *)nextdigit()) ) { // scan over the reduced list of digits 
643           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!!
644           switch (ineb ) {
645           case 0 :   // not a neighbour
646             break ;
647           case 1 :   // are neighbours 
648             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp(), digitN->GetId() ) ) ;
649             clusterdigitslist[iDigitInCluster] = digitN->GetIndexInList() ; 
650             iDigitInCluster++ ; 
651             digitsC->Remove(digitN) ;
652             break ;
653           case 2 :   // too far from each other
654             goto endofloop;   
655           } // switch
656           
657         } // while digitN
658         
659       endofloop: ;
660         nextdigit.Reset() ; 
661         
662       } // loop over cluster     
663
664     } // energy theshold  
665
666     
667   } // while digit
668
669   delete digitsC ;
670
671 }
672
673 //____________________________________________________________________________
674 void AliPHOSClusterizerv1::MakeUnfolding()
675 {
676   // Unfolds clusters using the shape of an ElectroMagnetic shower
677   // Performs unfolding of all EMC/CPV clusters
678
679   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
680   
681   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
682   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
683   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
684   TClonesArray * digits = gime->Digits() ; 
685   
686   // Unfold first EMC clusters 
687   if(fNumberOfEmcClusters > 0){
688
689     Int_t nModulesToUnfold = geom->GetNModules() ; 
690
691     Int_t numberofNotUnfolded = fNumberOfEmcClusters ; 
692     Int_t index ;   
693     for(index = 0 ; index < numberofNotUnfolded ; index++){
694       
695       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint *) emcRecPoints->At(index) ;
696       if(emcRecPoint->GetPHOSMod()> nModulesToUnfold)
697         break ;
698       
699       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
700       Int_t * maxAt = new Int_t[nMultipl] ;
701       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
702       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fEmcLocMaxCut,digits) ;
703       
704       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
705         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
706         emcRecPoints->Remove(emcRecPoint); 
707         emcRecPoints->Compress() ;
708         index-- ;
709         fNumberOfEmcClusters -- ;
710         numberofNotUnfolded-- ;
711       }
712       
713       delete[] maxAt ; 
714       delete[] maxAtEnergy ; 
715     }
716   } 
717   // Unfolding of EMC clusters finished
718
719
720   // Unfold now CPV clusters
721   if(fNumberOfCpvClusters > 0){
722     
723     Int_t nModulesToUnfold = geom->GetNModules() ;
724
725     Int_t numberofCpvNotUnfolded = fNumberOfCpvClusters ;     
726     Int_t index ;   
727     for(index = 0 ; index < numberofCpvNotUnfolded ; index++){
728       
729       AliPHOSRecPoint * recPoint = (AliPHOSRecPoint *) cpvRecPoints->At(index) ;
730
731       if(recPoint->GetPHOSMod()> nModulesToUnfold)
732         break ;
733       
734       AliPHOSEmcRecPoint * emcRecPoint = (AliPHOSEmcRecPoint*) recPoint ; 
735       
736       Int_t nMultipl = emcRecPoint->GetMultiplicity() ; 
737       Int_t * maxAt = new Int_t[nMultipl] ;
738       Float_t * maxAtEnergy = new Float_t[nMultipl] ;
739       Int_t nMax = emcRecPoint->GetNumberOfLocalMax(maxAt, maxAtEnergy,fCpvLocMaxCut,digits) ;
740       
741       if( nMax > 1 ) {     // if cluster is very flat (no pronounced maximum) then nMax = 0       
742         UnfoldCluster(emcRecPoint, nMax, maxAt, maxAtEnergy) ;
743         cpvRecPoints->Remove(emcRecPoint); 
744         cpvRecPoints->Compress() ;
745         index-- ;
746         numberofCpvNotUnfolded-- ;
747         fNumberOfCpvClusters-- ;
748       }
749       
750       delete[] maxAt ; 
751       delete[] maxAtEnergy ; 
752     } 
753   }
754   //Unfolding of Cpv clusters finished
755   
756 }
757
758 //____________________________________________________________________________
759 Double_t  AliPHOSClusterizerv1::ShowerShape(Double_t r)
760
761   // Shape of the shower (see PHOS TDR)
762   // If you change this function, change also the gradient evaluation in ChiSquare()
763
764   Double_t r4    = r*r*r*r ;
765   Double_t r295  = TMath::Power(r, 2.95) ;
766   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
767   return shape ;
768 }
769
770 //____________________________________________________________________________
771 void  AliPHOSClusterizerv1::UnfoldCluster(AliPHOSEmcRecPoint * iniEmc, 
772                                                  Int_t nMax, 
773                                                  int * maxAt, 
774                                                  Float_t * maxAtEnergy)
775
776   // Performs the unfolding of a cluster with nMax overlapping showers 
777
778   AliPHOSGetter * gime = AliPHOSGetter::GetInstance() ; 
779   const AliPHOSGeometry * geom = gime->PHOSGeometry() ;
780   const TClonesArray * digits = gime->Digits() ; 
781   TObjArray * emcRecPoints = gime->EmcRecPoints() ; 
782   TObjArray * cpvRecPoints = gime->CpvRecPoints() ; 
783
784   Int_t nPar = 3 * nMax ;
785   Float_t * fitparameters = new Float_t[nPar] ;
786
787   Bool_t rv = FindFit(iniEmc, maxAt, maxAtEnergy, nPar, fitparameters) ;
788   if( !rv ) {
789     // Fit failed, return and remove cluster
790     delete[] fitparameters ; 
791     return ;
792   }
793
794   // create ufolded rec points and fill them with new energy lists
795   // First calculate energy deposited in each sell in accordance with fit (without fluctuations): efit[]
796   // and later correct this number in acordance with actual energy deposition
797
798   Int_t nDigits = iniEmc->GetMultiplicity() ;  
799   Float_t * efit = new Float_t[nDigits] ;
800   Float_t xDigit=0.,zDigit=0.,distance=0. ;
801   Float_t xpar=0.,zpar=0.,epar=0.  ;
802   Int_t relid[4] ;
803   AliPHOSDigit * digit = 0 ;
804   Int_t * emcDigits = iniEmc->GetDigitsList() ;
805
806   Int_t iparam ;
807   Int_t iDigit ;
808   for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
809     digit = (AliPHOSDigit*) digits->At(emcDigits[iDigit] ) ;   
810     geom->AbsToRelNumbering(digit->GetId(), relid) ;
811     geom->RelPosInModule(relid, xDigit, zDigit) ;
812     efit[iDigit] = 0;
813
814     iparam = 0 ;    
815     while(iparam < nPar ){
816       xpar = fitparameters[iparam] ;
817       zpar = fitparameters[iparam+1] ;
818       epar = fitparameters[iparam+2] ;
819       iparam += 3 ;
820       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
821       distance =  TMath::Sqrt(distance) ;
822       efit[iDigit] += epar * ShowerShape(distance) ;
823     }
824   }
825   
826
827   // Now create new RecPoints and fill energy lists with efit corrected to fluctuations
828   // so that energy deposited in each cell is distributed betwin new clusters proportionally
829   // to its contribution to efit
830
831   Float_t * emcEnergies = iniEmc->GetEnergiesList() ;
832   Float_t ratio ;
833
834   iparam = 0 ;
835   while(iparam < nPar ){
836     xpar = fitparameters[iparam] ;
837     zpar = fitparameters[iparam+1] ;
838     epar = fitparameters[iparam+2] ;
839     iparam += 3 ;    
840     
841     AliPHOSEmcRecPoint * emcRP = 0 ;  
842
843     if(iniEmc->IsEmc()){ //create new entries in fEmcRecPoints...
844       
845       if(fNumberOfEmcClusters >= emcRecPoints->GetSize())
846         emcRecPoints->Expand(2*fNumberOfEmcClusters) ;
847       
848       (*emcRecPoints)[fNumberOfEmcClusters] = new AliPHOSEmcRecPoint("") ;
849       emcRP = (AliPHOSEmcRecPoint *) emcRecPoints->At(fNumberOfEmcClusters);
850       fNumberOfEmcClusters++ ;
851     }
852     else{//create new entries in fCpvRecPoints
853       if(fNumberOfCpvClusters >= cpvRecPoints->GetSize())
854         cpvRecPoints->Expand(2*fNumberOfCpvClusters) ;
855       
856       (*cpvRecPoints)[fNumberOfCpvClusters] = new AliPHOSCpvRecPoint("") ;
857       emcRP = (AliPHOSEmcRecPoint *) cpvRecPoints->At(fNumberOfCpvClusters);
858       fNumberOfCpvClusters++ ;
859     }
860     
861     Float_t eDigit ;
862     for(iDigit = 0 ; iDigit < nDigits ; iDigit ++){
863       digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ; 
864       geom->AbsToRelNumbering(digit->GetId(), relid) ;
865       geom->RelPosInModule(relid, xDigit, zDigit) ;
866       distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
867       distance =  TMath::Sqrt(distance) ;
868       ratio = epar * ShowerShape(distance) / efit[iDigit] ; 
869       eDigit = emcEnergies[iDigit] * ratio ;
870       emcRP->AddDigit( *digit, eDigit ) ;
871     }   
872   }
873  
874   delete[] fitparameters ; 
875   delete[] efit ; 
876   
877 }
878
879 //_____________________________________________________________________________
880 void AliPHOSClusterizerv1::UnfoldingChiSquare(Int_t & nPar, Double_t * Grad, Double_t & fret, Double_t * x, Int_t iflag)
881 {
882   // Calculates the Chi square for the cluster unfolding minimization
883   // Number of parameters, Gradient, Chi squared, parameters, what to do
884
885
886   TList * toMinuit = (TList*) gMinuit->GetObjectFit() ;
887
888   AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint*) toMinuit->At(0)  ;
889   TClonesArray * digits = (TClonesArray*)toMinuit->At(1)  ;
890
891
892   
893   //  AliPHOSEmcRecPoint * emcRP = (AliPHOSEmcRecPoint *) gMinuit->GetObjectFit() ; // EmcRecPoint to fit
894
895   Int_t * emcDigits     = emcRP->GetDigitsList() ;
896
897   Int_t nOdigits = emcRP->GetDigitsMultiplicity() ; 
898
899   Float_t * emcEnergies = emcRP->GetEnergiesList() ;
900
901   const AliPHOSGeometry * geom = AliPHOSGetter::GetInstance()->PHOSGeometry() ; 
902   fret = 0. ;     
903   Int_t iparam ;
904
905   if(iflag == 2)
906     for(iparam = 0 ; iparam < nPar ; iparam++)    
907       Grad[iparam] = 0 ; // Will evaluate gradient
908   
909   Double_t efit ;    
910
911   AliPHOSDigit * digit ;
912   Int_t iDigit ;
913
914   for( iDigit = 0 ; iDigit < nOdigits ; iDigit++) {
915
916     digit = (AliPHOSDigit*) digits->At( emcDigits[iDigit] ) ; 
917
918     Int_t relid[4] ;
919     Float_t xDigit ;
920     Float_t zDigit ;
921
922     geom->AbsToRelNumbering(digit->GetId(), relid) ;
923
924     geom->RelPosInModule(relid, xDigit, zDigit) ;
925
926      if(iflag == 2){  // calculate gradient
927        Int_t iParam = 0 ;
928        efit = 0 ;
929        while(iParam < nPar ){
930          Double_t distance = (xDigit - x[iParam]) * (xDigit - x[iParam]) ;
931          iParam++ ; 
932          distance += (zDigit - x[iParam]) * (zDigit - x[iParam]) ; 
933          distance = TMath::Sqrt( distance ) ; 
934          iParam++ ;      
935          efit += x[iParam] * ShowerShape(distance) ;
936          iParam++ ;
937        }
938        Double_t sum = 2. * (efit - emcEnergies[iDigit]) / emcEnergies[iDigit] ; // Here we assume, that sigma = sqrt(E) 
939        iParam = 0 ;
940        while(iParam < nPar ){
941          Double_t xpar = x[iParam] ;
942          Double_t zpar = x[iParam+1] ;
943          Double_t epar = x[iParam+2] ;
944          Double_t dr = TMath::Sqrt( (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar) );
945          Double_t shape = sum * ShowerShape(dr) ;
946          Double_t r4 = dr*dr*dr*dr ;
947          Double_t r295 = TMath::Power(dr,2.95) ;
948          Double_t deriv =-4. * dr*dr * ( 2.32 / ( (2.32 + 0.26 * r4) * (2.32 + 0.26 * r4) ) +
949                                          0.0316 * (1. + 0.0171 * r295) / ( ( 1. + 0.0652 * r295) * (1. + 0.0652 * r295) ) ) ;
950          
951          Grad[iParam] += epar * shape * deriv * (xpar - xDigit) ;  // Derivative over x    
952          iParam++ ; 
953          Grad[iParam] += epar * shape * deriv * (zpar - zDigit) ;  // Derivative over z         
954          iParam++ ; 
955          Grad[iParam] += shape ;                                  // Derivative over energy             
956          iParam++ ; 
957        }
958      }
959      efit = 0;
960      iparam = 0 ;
961
962      while(iparam < nPar ){
963        Double_t xpar = x[iparam] ;
964        Double_t zpar = x[iparam+1] ;
965        Double_t epar = x[iparam+2] ;
966        iparam += 3 ;
967        Double_t distance = (xDigit - xpar) * (xDigit - xpar) + (zDigit - zpar) * (zDigit - zpar)  ;
968        distance =  TMath::Sqrt(distance) ;
969        efit += epar * ShowerShape(distance) ;
970      }
971
972      fret += (efit-emcEnergies[iDigit])*(efit-emcEnergies[iDigit])/emcEnergies[iDigit] ; 
973      // Here we assume, that sigma = sqrt(E)
974   }
975
976 }
977
978 //____________________________________________________________________________
979 void AliPHOSClusterizerv1::Print(Option_t * option)const
980 {
981   // Print clusterizer parameters
982
983   if( strcmp(GetName(), "") !=0 ){
984     
985     // Print parameters
986  
987     TString taskName(GetName()) ; 
988     taskName.ReplaceAll(Version(), "") ;
989
990     cout << "---------------"<< taskName.Data() << " " << GetTitle()<< "-----------" << endl 
991          << "Clusterizing digits from the file: " << fHeaderFileName.Data() << endl 
992          << "                           Branch: " << fDigitsBranchTitle.Data() << endl 
993          << endl 
994          << "                       EMC Clustering threshold = " << fEmcClusteringThreshold << endl
995          << "                       EMC Local Maximum cut    = " << fEmcLocMaxCut << endl
996          << "                       EMC Logarothmic weight   = " << fW0 << endl
997          << endl
998          << "                       CPV Clustering threshold = " << fCpvClusteringThreshold << endl
999          << "                       CPV Local Maximum cut    = " << fCpvLocMaxCut << endl
1000        << "                       CPV Logarothmic weight   = " << fW0CPV << endl
1001          << endl ;
1002     if(fToUnfold)
1003       cout << " Unfolding on " << endl ;
1004     else
1005       cout << " Unfolding off " << endl ;
1006     
1007     cout << "------------------------------------------------------------------" <<endl ;
1008   }
1009   else
1010     cout << " AliPHOSClusterizerv1 not initialized " << endl ;
1011 }
1012 //____________________________________________________________________________
1013 void AliPHOSClusterizerv1::PrintRecPoints(Option_t * option)
1014 {
1015   // Prints list of RecPoints produced at the current pass of AliPHOSClusterizer
1016
1017   TObjArray * emcRecPoints = AliPHOSGetter::GetInstance()->EmcRecPoints() ; 
1018   TObjArray * cpvRecPoints = AliPHOSGetter::GetInstance()->CpvRecPoints() ; 
1019
1020   cout << "AliPHOSClusterizerv1: : event "<<gAlice->GetEvNumber() << endl ;
1021   cout << "       Found "<< emcRecPoints->GetEntriesFast() << " EMC Rec Points and " 
1022            << cpvRecPoints->GetEntriesFast() << " CPV RecPoints" << endl ;
1023
1024   fRecPointsInRun +=  emcRecPoints->GetEntriesFast() ; 
1025   fRecPointsInRun +=  cpvRecPoints->GetEntriesFast() ; 
1026
1027   if(strstr(option,"all")) {
1028     cout << "EMC clusters " << endl ;
1029     cout << " Index  Ene(MeV)   Multi  Module     X      Y      Z    Lambda 1   Lambda 2  # of prim  Primaries list "      <<  endl;      
1030     
1031     Int_t index ;
1032     for (index = 0 ; index < emcRecPoints->GetEntries() ; index++) {
1033       AliPHOSEmcRecPoint * rp = (AliPHOSEmcRecPoint * )emcRecPoints->At(index) ; 
1034       TVector3  locpos;  
1035       rp->GetLocalPosition(locpos);
1036       Float_t lambda[2]; 
1037       rp->GetElipsAxis(lambda);
1038       Int_t * primaries; 
1039       Int_t nprimaries;
1040       primaries = rp->GetPrimaries(nprimaries);
1041
1042       cout << setw(4) << rp->GetIndexInList() << "   " 
1043            << setw(7) << setprecision(3) << rp->GetEnergy() << "           " 
1044            << setw(3) <<         rp->GetMultiplicity() << "  " 
1045            << setw(1) <<              rp->GetPHOSMod() << "  " 
1046            << setw(6) << setprecision(2) << locpos.X() << "  " 
1047            << setw(6) << setprecision(2) << locpos.Y() << "  " 
1048            << setw(6) << setprecision(2) << locpos.Z() << "  "
1049            << setw(4) << setprecision(2) << lambda[0]  << "  "
1050            << setw(4) << setprecision(2) << lambda[1]  << "  "
1051            << setw(2) << nprimaries << "  " ;
1052      
1053       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1054         cout << setw(4) <<   primaries[iprimary] << "  "  ;
1055       cout << endl ;     
1056     }
1057
1058     //Now plot CPV recPoints
1059     cout << "EMC clusters " << endl ;
1060     cout << "  Index    " 
1061          << "  Multi    "
1062          << "  Module   "  
1063          << "    X      "
1064          << "    Y      "
1065          << "    Z      "
1066          << " # of prim "
1067          << " Primaries list "      <<  endl;      
1068     
1069     for (index = 0 ; index < cpvRecPoints->GetEntries() ; index++) {
1070       AliPHOSRecPoint * rp = (AliPHOSRecPoint * )cpvRecPoints->At(index) ; 
1071       cout << setw(6) << rp->GetIndexInList() << "     ";
1072       cout << setw(6) << rp->GetPHOSMod()     << "        CPV     ";
1073       
1074       TVector3  locpos;  
1075       rp->GetLocalPosition(locpos);
1076       cout << setw(6) <<  locpos.X()          << "     ";
1077       cout << setw(6) <<  locpos.Y()          << "     ";
1078       cout << setw(6) <<  locpos.Z()          << "     ";
1079       
1080       Int_t * primaries; 
1081       Int_t nprimaries ; 
1082       primaries = rp->GetPrimaries(nprimaries);
1083       cout << setw(6) <<    nprimaries         << "     ";
1084
1085       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++)
1086         cout << setw(4)  <<  primaries[iprimary] << " ";
1087       cout << endl;      
1088     }
1089
1090
1091     cout << "-----------------------------------------------------------------------"<<endl ;
1092   }
1093 }
1094