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