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