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