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