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