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