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