Including assert.h (Solaris)
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALClusterizerv1.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 /* $Id$ */
18
19 //*-- Author: Yves Schutz (SUBATECH)  & Dmitri Peressounko (SUBATECH & Kurchatov Institute)
20 //  August 2002 Yves Schutz: clone PHOS as closely as possible and intoduction
21 //                           of new  IO (à la PHOS)
22 //////////////////////////////////////////////////////////////////////////////
23 //  Clusterization class. Performs clusterization (collects neighbouring active cells) and 
24 //  unfolds the clusters having several local maxima.  
25 //  Results are stored in TreeR#, branches EMCALTowerRP (EMC recPoints),
26 //  EMCALPreShoRP (CPV RecPoints) and AliEMCALClusterizer (Clusterizer with all 
27 //  parameters including input digits branch title, thresholds etc.)
28 //  This TTask is normally called from Reconstructioner, but can as well be used in 
29 //  standalone mode.
30 // Use Case:
31 //  root [0] AliEMCALClusterizerv1 * cl = new AliEMCALClusterizerv1("galice.root")  
32 //  Warning in <TDatabasePDG::TDatabasePDG>: object already instantiated
33 //               //reads gAlice from header file "..."                      
34 //  root [1] cl->ExecuteTask()  
35 //               //finds RecPoints in all events stored in galice.root
36 //  root [2] cl->SetDigitsBranch("digits2") 
37 //               //sets another title for Digitis (input) branch
38 //  root [3] cl->SetRecPointsBranch("recp2")  
39 //               //sets another title four output branches
40 //  root [4] cl->SetTowerLocalMaxCut(0.03)  
41 //               //set clusterization parameters
42 //  root [5] cl->ExecuteTask("deb all time")  
43 //               //once more finds RecPoints options are 
44 //               // deb - print number of found rec points
45 //               // deb all - print number of found RecPoints and some their characteristics 
46 //               // time - print benchmarking results
47
48 // --- ROOT system ---
49
50 #include "TROOT.h" 
51 #include "TFile.h" 
52 #include "TFolder.h" 
53 #include "TMath.h" 
54 #include "TMinuit.h"
55 #include "TTree.h" 
56 #include "TSystem.h" 
57 #include "TBenchmark.h"
58
59 // --- Standard library ---
60
61
62 // --- AliRoot header files ---
63 #include "AliRunLoader.h"
64 #include "AliRun.h"
65 #include "AliEMCALLoader.h"
66 #include "AliEMCALClusterizerv1.h"
67 #include "AliEMCALRecPoint.h"
68 #include "AliEMCALDigit.h"
69 #include "AliEMCALDigitizer.h"
70 #include "AliEMCAL.h"
71 #include "AliEMCALGeometry.h"
72
73 ClassImp(AliEMCALClusterizerv1)
74
75 Int_t addOn[20][60][60]; 
76
77 //____________________________________________________________________________
78   AliEMCALClusterizerv1::AliEMCALClusterizerv1() : AliEMCALClusterizer()
79 {
80   // default ctor (to be used mainly by Streamer)
81   
82   InitParameters() ; 
83   fDefaultInit = kTRUE ;
84   for(Int_t is=0;is<20;is++){ 
85     for(Int_t i=0;i<60;i++){ 
86       for(Int_t j=0;j<60;j++){ 
87         addOn[is][i][j]=0;
88       }
89     }
90   }
91 //PH   cout<<"file to read 1"<<endl;
92   ReadFile();
93 //PH   cout<<"file read 1"<<endl;
94 }
95
96 //____________________________________________________________________________
97 AliEMCALClusterizerv1::AliEMCALClusterizerv1(const TString alirunFileName, const TString eventFolderName)
98 :AliEMCALClusterizer(alirunFileName, eventFolderName)
99 {
100   // ctor with the indication of the file where header Tree and digits Tree are stored
101   
102   InitParameters() ; 
103   Init() ;
104   fDefaultInit = kFALSE ; 
105   for(Int_t is=0;is<20;is++){ 
106     for(Int_t i=0;i<60;i++){ 
107       for(Int_t j=0;j<60;j++){ 
108         addOn[is][i][j]=0;
109       }
110     }
111   }
112 //PH   cout<<"file to read 2"<<endl;
113   ReadFile();
114 //PH   cout<<"file read 2"<<endl;
115
116 }
117
118 //____________________________________________________________________________
119   AliEMCALClusterizerv1::~AliEMCALClusterizerv1()
120 {
121   // dtor
122   
123 }
124
125 //____________________________________________________________________________
126 const TString AliEMCALClusterizerv1::BranchName() const 
127
128    return GetName();
129
130 }
131
132 //____________________________________________________________________________
133 Float_t  AliEMCALClusterizerv1::Calibrate(Int_t amp) const
134 {
135   //To be replased later by the method, reading individual parameters from the database 
136   return -fADCpedestalECA + amp * fADCchannelECA ; 
137 }
138
139 //____________________________________________________________________________
140 void AliEMCALClusterizerv1::Exec(Option_t * option)
141 {
142   // Steering method to perform clusterization for events
143   // in the range from fFirstEvent to fLastEvent.
144   // This range is optionally set by SetEventRange().
145   // if fLastEvent=-1 (by default), then process events until the end.
146
147   if(strstr(option,"tim"))
148     gBenchmark->Start("EMCALClusterizer"); 
149   
150   if(strstr(option,"print"))
151     Print("") ; 
152
153   //AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
154   AliRunLoader *rl = AliRunLoader::GetRunLoader();
155   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(rl->GetDetectorLoader("EMCAL"));
156
157   if (fLastEvent == -1) 
158     fLastEvent = rl->GetNumberOfEvents() - 1;
159   Int_t nEvents   = fLastEvent - fFirstEvent + 1;
160
161   Int_t ievent ;
162   rl->LoadDigits("EMCAL");
163   for (ievent = fFirstEvent; ievent <= fLastEvent; ievent++) {
164     rl->GetEvent(ievent);
165     GetCalibrationParameters() ;
166
167     fNumberOfECAClusters = 0;
168            
169     MakeClusters() ;
170
171     if(fToUnfold)
172       MakeUnfolding() ;
173
174     WriteRecPoints() ;
175
176     if(strstr(option,"deb"))  
177       PrintRecPoints(option) ;
178
179     //increment the total number of recpoints per run   
180     fRecPointsInRun += emcalLoader->RecPoints()->GetEntriesFast() ;  
181   }
182   
183   Unload();
184
185   if(strstr(option,"tim")){
186     gBenchmark->Stop("EMCALClusterizer");
187     printf("Exec took %f seconds for Clusterizing %f seconds per event", 
188          gBenchmark->GetCpuTime("EMCALClusterizer"), gBenchmark->GetCpuTime("EMCALClusterizer")/nEvents ) ;
189   }
190
191   SaveHists();
192
193 }
194
195 //____________________________________________________________________________
196 Bool_t AliEMCALClusterizerv1::FindFit(AliEMCALRecPoint * emcRP, AliEMCALDigit ** maxAt, Float_t * maxAtEnergy,
197                                     Int_t nPar, Float_t * fitparameters) const
198
199   // Calls TMinuit to fit the energy distribution of a cluster with several maxima 
200   // The initial values for fitting procedure are set equal to the positions of local maxima.
201   // Cluster will be fitted as a superposition of nPar/3 electromagnetic showers
202
203   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
204   TClonesArray *digits = emcalLoader->Digits();
205   /*
206   AliEMCALGetter * gime = AliEMCALGetter::Instance() ; 
207   TClonesArray * digits = gime->Digits() ; 
208   */
209
210   gMinuit->mncler();                     // Reset Minuit's list of paramters
211   gMinuit->SetPrintLevel(-1) ;           // No Printout
212   gMinuit->SetFCN(AliEMCALClusterizerv1::UnfoldingChiSquare) ;  
213                                          // To set the address of the minimization function 
214   TList * toMinuit = new TList();
215   toMinuit->AddAt(emcRP,0) ;
216   toMinuit->AddAt(digits,1) ;
217   
218   gMinuit->SetObjectFit(toMinuit) ;         // To tranfer pointer to UnfoldingChiSquare
219
220   // filling initial values for fit parameters
221   AliEMCALDigit * digit ;
222
223   Int_t ierflg  = 0; 
224   Int_t index   = 0 ;
225   Int_t nDigits = (Int_t) nPar / 3 ;
226
227   Int_t iDigit ;
228
229   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance(); //gime->EMCALGeometry() ; 
230
231   for(iDigit = 0; iDigit < nDigits; iDigit++){
232     digit = maxAt[iDigit]; 
233
234     Int_t relid[2] ;
235     Float_t x = 0.;
236     Float_t z = 0.;
237     geom->AbsToRelNumbering(digit->GetId(), relid) ;
238     geom->PosInAlice(relid, x, z) ;
239
240     Float_t energy = maxAtEnergy[iDigit] ;
241
242     gMinuit->mnparm(index, "x",  x, 0.1, 0, 0, ierflg) ;
243     index++ ;   
244     if(ierflg != 0){ 
245       Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : x = %f",  x ) ;
246       return kFALSE;
247     }
248     gMinuit->mnparm(index, "z",  z, 0.1, 0, 0, ierflg) ;
249     index++ ;   
250     if(ierflg != 0){
251        Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : z = %f", z) ;
252       return kFALSE;
253     }
254     gMinuit->mnparm(index, "Energy",  energy , 0.05*energy, 0., 4.*energy, ierflg) ;
255     index++ ;   
256     if(ierflg != 0){
257      Error("FindFit", "EMCAL Unfolding  Unable to set initial value for fit procedure : energy = %f", energy) ;      
258       return kFALSE;
259     }
260   }
261
262   Double_t p0 = 0.1 ; // "Tolerance" Evaluation stops when EDM = 0.0001*p0 ; The number of function call slightly
263                       //  depends on it. 
264   Double_t p1 = 1.0 ;
265   Double_t p2 = 0.0 ;
266
267   gMinuit->mnexcm("SET STR", &p2, 0, ierflg) ;   // force TMinuit to reduce function calls  
268   gMinuit->mnexcm("SET GRA", &p1, 1, ierflg) ;   // force TMinuit to use my gradient  
269   gMinuit->SetMaxIterations(5);
270   gMinuit->mnexcm("SET NOW", &p2 , 0, ierflg) ;  // No Warnings
271
272   gMinuit->mnexcm("MIGRAD", &p0, 0, ierflg) ;    // minimize 
273
274   if(ierflg == 4){  // Minimum not found   
275     Error("FindFit", "EMCAL Unfolding  Fit not converged, cluster abandoned " ) ;      
276     return kFALSE ;
277   }            
278   for(index = 0; index < nPar; index++){
279     Double_t err ;
280     Double_t val ;
281     gMinuit->GetParameter(index, val, err) ;    // Returns value and error of parameter index
282     fitparameters[index] = val ;
283    }
284
285   delete toMinuit ;
286   return kTRUE;
287
288 }
289
290 //____________________________________________________________________________
291 void AliEMCALClusterizerv1::GetCalibrationParameters() 
292 {
293   // Gets the parameters for the calibration from the digitizer
294   //AliEMCALGetter * gime = AliEMCALGetter::Instance() ;
295   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
296
297   if ( !emcalLoader->Digitizer() ) 
298     emcalLoader->LoadDigitizer();
299   AliEMCALDigitizer * dig = dynamic_cast<AliEMCALDigitizer*>(emcalLoader->Digitizer());
300
301   fADCchannelECA   = dig->GetECAchannel() ;
302   fADCpedestalECA  = dig->GetECApedestal();
303 //PH  cout<<"ChannelECA, peds "<<fADCchannelECA<<" "<<fADCpedestalECA<<endl;
304 }
305
306 //____________________________________________________________________________
307 void AliEMCALClusterizerv1::Init()
308 {
309   // Make all memory allocations which can not be done in default constructor.
310   // Attach the Clusterizer task to the list of EMCAL tasks
311   
312   //AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance() ;
313   AliRunLoader *rl = AliRunLoader::GetRunLoader();
314   AliEMCALGeometry *geom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
315   AliInfo(Form("geom 0x%x",geom));
316
317 //Sub  fNTowers = geom->GetNZ() *  geom->GetNPhi() ;
318   fNTowers =400;
319   if(!gMinuit) 
320     gMinuit = new TMinuit(100) ;
321  //Sub if ( !gime->Clusterizer() ) 
322  //Sub   gime->PostClusterizer(this); 
323  BookHists();
324 }
325
326 //____________________________________________________________________________
327 void AliEMCALClusterizerv1::InitParameters()
328
329   // Initializes the parameters for the Clusterizer
330   fNumberOfECAClusters = 0;
331
332   fECAClusteringThreshold   = 0.5;  // value obtained from Alexei
333   fECALocMaxCut = 0.03;
334
335   fECAW0     = 4.5 ;
336   fTimeGate = 1.e-8 ; 
337   fToUnfold = kFALSE ;
338   fRecPointsInRun  = 0 ;
339   fMinECut = 0.3;
340 }
341
342 //____________________________________________________________________________
343 Int_t AliEMCALClusterizerv1::AreNeighbours(AliEMCALDigit * d1, AliEMCALDigit * d2)const
344 {
345   // Gives the neighbourness of two digits = 0 are not neighbour but continue searching 
346   //                                       = 1 are neighbour
347   //                                       = 2 are not neighbour but do not continue searching
348   // neighbours are defined as digits having at least a common vertex 
349   // The order of d1 and d2 is important: first (d1) should be a digit already in a cluster 
350   //                                      which is compared to a digit (d2)  not yet in a cluster  
351
352   //AliEMCALGeometry * geom = AliEMCALGetter::Instance()->EMCALGeometry() ;
353    AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance() ;
354
355   Int_t rv = 0 ; 
356
357   Int_t relid1[2] ; 
358  //Sub geom->AbsToRelNumbering(d1->GetId(), relid1) ; 
359     Int_t nSupMod=0;
360     Int_t nTower=0;
361     Int_t nIphi=0;
362     Int_t nIeta=0;
363     Int_t iphi=0;
364     Int_t ieta=0;
365    geom->GetCellIndex(d1->GetId(), nSupMod,nTower,nIphi,nIeta);
366    geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
367    relid1[0]=ieta;
368    relid1[1]=iphi;
369
370
371     Int_t nSupMod1=0;
372     Int_t nTower1=0;
373     Int_t nIphi1=0;
374     Int_t nIeta1=0;
375     Int_t iphi1=0;
376     Int_t ieta1=0;
377    geom->GetCellIndex(d2->GetId(), nSupMod1,nTower1,nIphi1,nIeta1);
378    geom->GetCellPhiEtaIndexInSModule(nSupMod1, nTower1,nIphi1,nIeta1, iphi1,ieta1);
379    Int_t relid2[2] ; 
380    relid2[0]=ieta1;
381    relid2[1]=iphi1;
382
383   //Sub  geom->AbsToRelNumbering(d2->GetId(), relid2) ; 
384   
385
386   Int_t rowdiff = TMath::Abs( relid1[0] - relid2[0] ) ;  
387   Int_t coldiff = TMath::Abs( relid1[1] - relid2[1] ) ;  
388   
389   if (( coldiff <= 1 )  && ( rowdiff <= 1 )){
390       rv = 1 ; 
391   }
392   else {
393     if((relid2[0] > relid1[0]) && (relid2[1] > relid1[1]+1)) 
394       rv = 2; //  Difference in row numbers is too large to look further 
395   }
396  
397   if (gDebug == 2 ) 
398 if(rv==1)printf("AreNeighbours: neighbours=%d, id1=%d, relid1=%d,%d \n id2=%d, relid2=%d,%d \n", 
399          rv, d1->GetId(), relid1[0], relid1[1],
400          d2->GetId(), relid2[0], relid2[1]) ;   
401   
402   return rv ; 
403 }
404
405 //____________________________________________________________________________
406 void AliEMCALClusterizerv1::Unload() 
407 {
408   // Unloads the Digits and RecPoints
409   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
410     
411   emcalLoader->UnloadDigits() ; 
412   emcalLoader->UnloadRecPoints() ; 
413 }
414  
415 //____________________________________________________________________________
416 void AliEMCALClusterizerv1::WriteRecPoints()
417 {
418
419   // Creates new branches with given title
420   // fills and writes into TreeR.
421
422   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
423
424   TObjArray * aECARecPoints = emcalLoader->RecPoints() ; 
425
426   TClonesArray * digits = emcalLoader->Digits() ; 
427   TTree * treeR = emcalLoader->TreeR();  
428   if ( treeR==0 ) {
429     emcalLoader->MakeRecPointsContainer();
430     treeR = emcalLoader->TreeR();  
431   }
432   Int_t index ;
433
434   //Evaluate position, dispersion and other RecPoint properties for EC section
435   for(index = 0; index < aECARecPoints->GetEntries(); index++)
436     (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->EvalAll(fECAW0,digits) ;
437   
438   aECARecPoints->Sort() ;
439
440   for(index = 0; index < aECARecPoints->GetEntries(); index++)
441     (dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(index)))->SetIndexInList(index) ;
442
443   Int_t bufferSize = 32000 ;    
444   Int_t splitlevel = 0 ; 
445
446   //EC section branch
447   
448   TBranch * branchECA = 0;
449   if ((branchECA = treeR->GetBranch("EMCALECARP")))
450     branchECA->SetAddress(&aECARecPoints);
451   else
452     treeR->Branch("EMCALECARP","TObjArray",&aECARecPoints,bufferSize,splitlevel);
453   treeR->Fill() ;
454
455   emcalLoader->WriteRecPoints("OVERWRITE");
456   //emcalLoader->WriteClusterizer("OVERWRITE");
457   emcalLoader->WriteReconstructioner("OVERWRITE");
458 }
459
460 //____________________________________________________________________________
461 void AliEMCALClusterizerv1::MakeClusters()
462 {
463   // Steering method to construct the clusters stored in a list of Reconstructed Points
464   // A cluster is defined as a list of neighbour digits
465
466   cout << "Star of " << __PRETTY_FUNCTION__ << endl;
467   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
468
469   TObjArray * aECARecPoints = emcalLoader->RecPoints() ; 
470     
471   AliEMCALGeometry * geom = AliEMCALGeometry::GetInstance();
472   if (geom==0) 
473      AliFatal("Did not get geometry from EMCALLoader");
474
475   aECARecPoints->Clear();
476
477   TClonesArray * digits = emcalLoader->Digits() ; 
478   TClonesArray * digitsC =  dynamic_cast<TClonesArray*>(digits->Clone()) ;
479
480   // Clusterization starts    
481   TIter nextdigit(digitsC) ;
482   AliEMCALDigit * digit;
483
484   cout << "Outer Loop" << endl;
485   while ( (digit = dynamic_cast<AliEMCALDigit *>(nextdigit())) ) { // scan over the list of digitsC
486     AliEMCALRecPoint * clu = 0 ; 
487     
488     TArrayI clusterECAdigitslist(50000);   
489
490 ////////////////////////temp solution, embedding///////////////////
491    int nSupMod=0, nTower=0, nIphi=0, nIeta=0;
492    int iphi=0, ieta=0;
493    geom->GetCellIndex(digit->GetId(), nSupMod,nTower,nIphi,nIeta);
494    geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
495    
496    //cout << "Some checking" << endl;
497    if (geom->CheckAbsCellId(digit->GetId()) && (Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) > fECAClusteringThreshold  ) ){
498       //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"11 digit, add "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
499 //          cout<<"crossed the threshold "<<endl;
500       Int_t iDigitInECACluster = 0;
501       // start a new Tower RecPoint
502       if(fNumberOfECAClusters >= aECARecPoints->GetSize()) 
503           aECARecPoints->Expand(2*fNumberOfECAClusters+1) ;
504       AliEMCALRecPoint * rp = new  AliEMCALRecPoint("") ; 
505       //rp->SetGeom(geom);
506       aECARecPoints->AddAt(rp, fNumberOfECAClusters) ;
507       clu = dynamic_cast<AliEMCALRecPoint *>(aECARecPoints->At(fNumberOfECAClusters)) ; 
508       fNumberOfECAClusters++ ; 
509       clu->AddDigit(*digit, Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1])) ; 
510       clusterECAdigitslist[iDigitInECACluster] = digit->GetIndexInList() ;      
511       iDigitInECACluster++ ; 
512 //    cout<<" 1st setp:cluno, digNo "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
513       digitsC->Remove(digit) ; 
514       AliDebug(1,Form("MakeClusters: OK id = %d, ene = %f , thre = %f \n", digit->GetId(),Calibrate(digit->GetAmp()), fECAClusteringThreshold));  
515       nextdigit.Reset() ;
516       
517       AliEMCALDigit * digitN ; 
518       Int_t index = 0 ;
519
520       // Find the neighbours
521       while (index < iDigitInECACluster){ // scan over digits already in cluster 
522         digit =  (AliEMCALDigit*)digits->At(clusterECAdigitslist[index]);
523         index++ ; 
524         while ( (digitN = (AliEMCALDigit *)nextdigit())) { // scan over the reduced list of digits 
525           // cout<<"we have new digit "<<endl;
526           // check that the digit is above the min E Cut
527           geom->GetCellIndex(digitN->GetId(), nSupMod,nTower,nIphi,nIeta);
528           geom->GetCellPhiEtaIndexInSModule(nSupMod,nTower,nIphi,nIeta, iphi,ieta);
529           if( Calibrate(digitN->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) < fMinECut  )  digitsC->Remove(digitN);
530           //cout<<" new digit above ECut "<<endl;
531           Int_t ineb = AreNeighbours(digit, digitN);       // call (digit,digitN) in THAT oder !!!!! 
532 //      cout<<" new digit neighbour?? "<<ineb<<endl;
533          switch (ineb ) {
534           case 0 :   // not a neighbour
535             break ;
536           case 1 :   // are neighbours 
537             //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"22 digit, add "<<nSupMod<<" "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
538             clu->AddDigit(*digitN, Calibrate( digitN->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) ) ;
539             clusterECAdigitslist[iDigitInECACluster] = digitN->GetIndexInList() ; 
540             iDigitInECACluster++ ; 
541 //    cout<<"when neighbour: cluno, digNo "<<digit->GetId()<<" "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
542             digitsC->Remove(digitN) ;
543 //          break ;
544 //          case 2 :   // too far from each other
545 //Subh      goto endofloop1;   
546 //          cout<<"earlier go to end of loop"<<endl;   
547           } // switch
548     //cout<<"in nextDigit loop "<<fNumberOfECAClusters<<" "<<iDigitInECACluster<<endl;
549         } // while digitN
550         //cout << "Done neighbout searching" << endl;
551         nextdigit.Reset() ; 
552       } // loop over ECA cluster
553     } // energy threshold
554     else if(Calibrate(digit->GetAmp()+addOn[nSupMod-1][ieta-1][iphi-1]) < fMinECut  ){
555       //if(addOn[nSupMod-1][ieta-1][iphi-1]>0)cout<<"33 digit, add  "<<ieta<<" "<<iphi<<" "<<addOn[nSupMod-1][ieta-1][iphi-1]<<" "<<digit->GetAmp()<<endl;
556       digitsC->Remove(digit);
557     }
558     //cout<<"after endofloop: cluno, digNo "<<fNumberOfECAClusters<<endl;
559   } // while digit  
560   delete digitsC ;
561   AliDebug(1,Form("total no of clusters %d from %d digits",fNumberOfECAClusters,digits->GetEntriesFast())); 
562 }
563
564 //____________________________________________________________________________
565 void AliEMCALClusterizerv1::MakeUnfolding() const
566 {
567   Fatal("AliEMCALClusterizerv1::MakeUnfolding", "--> Unfolding not implemented") ;
568  
569 }
570
571 //____________________________________________________________________________
572 Double_t  AliEMCALClusterizerv1::ShowerShape(Double_t r)
573
574   // Shape of the shower (see EMCAL TDR)
575   // If you change this function, change also the gradient evaluation in ChiSquare()
576
577   Double_t r4    = r*r*r*r ;
578   Double_t r295  = TMath::Power(r, 2.95) ;
579   Double_t shape = TMath::Exp( -r4 * (1. / (2.32 + 0.26 * r4) + 0.0316 / (1 + 0.0652 * r295) ) ) ;
580   return shape ;
581 }
582
583 //____________________________________________________________________________
584 void  AliEMCALClusterizerv1::UnfoldCluster(AliEMCALRecPoint * /*iniTower*/, 
585                                            Int_t /*nMax*/, 
586                                            AliEMCALDigit ** /*maxAt*/, 
587                                            Float_t * /*maxAtEnergy*/) const
588 {
589   // Performs the unfolding of a cluster with nMax overlapping showers 
590   
591   Fatal("UnfoldCluster", "--> Unfolding not implemented") ;
592
593 }
594
595 //_____________________________________________________________________________
596 void AliEMCALClusterizerv1::UnfoldingChiSquare(Int_t & /*nPar*/, Double_t * /*Grad*/,
597                                                Double_t & /*fret*/,
598                                                Double_t * /*x*/, Int_t /*iflag*/)
599 {
600   // Calculates the Chi square for the cluster unfolding minimization
601   // Number of parameters, Gradient, Chi squared, parameters, what to do
602   
603   ::Fatal("UnfoldingChiSquare","Unfolding not implemented") ;
604 }
605 //____________________________________________________________________________
606 void AliEMCALClusterizerv1::Print(Option_t * /*option*/)const
607 {
608   // Print clusterizer parameters
609
610   TString message("\n") ; 
611   
612   if( strcmp(GetName(), "") !=0 ){
613     
614     // Print parameters
615  
616     TString taskName(GetName()) ;
617     taskName.ReplaceAll(Version(), "") ;
618     
619     printf("--------------- "); 
620     printf(taskName.Data()) ; 
621     printf(" "); 
622     printf(GetTitle()) ; 
623     printf("-----------\n");  
624     printf("Clusterizing digits from the file: "); 
625     printf(taskName.Data());  
626     printf("\n                           Branch: "); 
627     printf(GetName()); 
628     printf("\n                       ECA Local Maximum cut    = %f", fECALocMaxCut); 
629     printf("\n                       ECA Logarithmic weight   = %f", fECAW0); 
630     if(fToUnfold)
631       printf("\nUnfolding on\n");
632     else
633       printf("\nUnfolding off\n");
634     
635     printf("------------------------------------------------------------------"); 
636   }
637   else
638     printf("AliEMCALClusterizerv1 not initialized ") ;
639 }
640
641 //____________________________________________________________________________
642 void AliEMCALClusterizerv1::PrintRecPoints(Option_t * option)
643 {
644   // Prints list of RecPoints produced at the current pass of AliEMCALClusterizer
645   AliEMCALLoader *emcalLoader = dynamic_cast<AliEMCALLoader*>(AliRunLoader::GetRunLoader()->GetDetectorLoader("EMCAL"));
646   TObjArray * aECARecPoints = emcalLoader->RecPoints() ; 
647     
648   printf("PrintRecPoints: Clusterization result:") ; 
649   
650   printf("event # %d\n", emcalLoader->GetRunLoader()->GetEventNumber() ) ;
651   printf("           Found %d ECA Rec Points\n ", 
652          aECARecPoints->GetEntriesFast()) ; 
653
654   fRecPointsInRun +=  aECARecPoints->GetEntriesFast() ; 
655   
656   if(strstr(option,"all")) {
657     Int_t index =0;
658     printf("\n-----------------------------------------------------------------------\n") ;
659     printf("Clusters in ECAL section\n") ;
660     printf("Index    Ene(GeV) Multi Module     phi     r   theta    X    Y      Z   Dispersion Lambda 1   Lambda 2  # of prim  Primaries list\n") ;      
661    Float_t maxE=0; 
662    Float_t maxL1=0; 
663    Float_t maxL2=0; 
664    Float_t maxDis=0; 
665
666     for (index = 0 ; index < aECARecPoints->GetEntries() ; index++) {
667       AliEMCALRecPoint * rp = dynamic_cast<AliEMCALRecPoint * >(aECARecPoints->At(index)) ; 
668       TVector3  globalpos;  
669       //rp->GetGlobalPosition(globalpos);
670       TVector3  localpos;  
671       rp->GetLocalPosition(localpos);
672       Float_t lambda[2]; 
673       rp->GetElipsAxis(lambda);
674       Int_t * primaries; 
675       Int_t nprimaries;
676       primaries = rp->GetPrimaries(nprimaries);
677       printf("\n%6d  %8.4f  %3d     %4.1f    %4.1f %4.1f  %4.1f %4.1f %4.1f    %4.1f   %4f  %4f    %2d     : ", 
678              rp->GetIndexInList(), rp->GetEnergy(), rp->GetMultiplicity(),
679              globalpos.X(), globalpos.Y(), globalpos.Z(), localpos.X(), localpos.Y(), localpos.Z(), 
680              rp->GetDispersion(), lambda[0], lambda[1], nprimaries) ; 
681   /////////////
682       if(rp->GetEnergy()>maxE){
683               maxE=rp->GetEnergy();
684               maxL1=lambda[0];
685               maxL2=lambda[1];
686               maxDis=rp->GetDispersion();
687       }
688       pointE->Fill(rp->GetEnergy());
689       pointL1->Fill(lambda[0]);
690       pointL2->Fill(lambda[1]);
691       pointDis->Fill(rp->GetDispersion());
692       pointMult->Fill(rp->GetMultiplicity());
693       ///////////// 
694       for (Int_t iprimary=0; iprimary<nprimaries; iprimary++) {
695         printf("%d ", primaries[iprimary] ) ; 
696       } 
697     }
698
699       MaxE->Fill(maxE);
700       MaxL1->Fill(maxL1);
701       MaxL2->Fill(maxL2);
702       MaxDis->Fill(maxDis);
703
704
705     printf("\n-----------------------------------------------------------------------\n");
706   }
707 }
708   void AliEMCALClusterizerv1::BookHists()
709 {
710         pointE = new TH1F("pointE","point energy", 2000, 0.0, 150.);
711         pointL1 = new TH1F("pointL1","point L1", 1000, 0.0, 3.);
712         pointL2 = new TH1F("pointL2","point L2", 1000, 0.0, 3.);
713         pointDis = new TH1F("pointDis","point Dis", 1000, 0.0, 3.);
714         pointMult = new TH1F("pointMult","point Mult", 100, 0.0, 100.);
715         digitAmp = new TH1F("digitAmp","Digit Amplitude", 2000, 0.0, 5000.);
716         MaxE = new TH1F("maxE","Max point energy", 2000, 0.0, 150.);
717         MaxL1 = new TH1F("maxL1","Max point L1", 1000, 0.0, 3.);
718         MaxL2 = new TH1F("maxL2","Max point L2", 1000, 0.0, 3.);
719         MaxDis = new TH1F("maxDis","Max point Dis", 1000, 0.0, 3.);
720 }
721 void AliEMCALClusterizerv1::SaveHists()
722 {
723  recofile=new TFile("reco.root","RECREATE"); 
724         pointE->Write();
725         pointL1->Write();
726         pointL2->Write();
727         pointDis->Write();
728         pointMult->Write();
729         digitAmp->Write();
730       MaxE->Write();
731       MaxL1->Write();
732       MaxL2->Write();
733       MaxDis->Write();
734 recofile->Close();
735 }
736
737 void AliEMCALClusterizerv1::ReadFile()
738 {
739   return; // 3-jan-05
740         FILE *fp = fopen("hijing1.dat","r");
741         for(Int_t line=0;line<9113;line++){
742          Int_t eg,l1,l2,sm;
743          Int_t ncols0;
744          ncols0 = fscanf(fp,"%d %d %d %d",&sm,&l1,&l2,&eg);
745         // cout<<eg<<" "<<l1<<" "<<l2<<endl;
746          addOn[sm-1][l1-1][l2-1]=eg;
747          //addOn[sm-1][l1-1][l2-1]=0;
748         }
749 }
750