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