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