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