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