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