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