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