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