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