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