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