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