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