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