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