]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
bad channels up to LHC12d added, Marcel
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALReconstructor.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 //_________________________________________________________________________
19 //--
20 //-- Yves Schutz (SUBATECH) 
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and 
22 // derived from STEER/AliReconstructor. 
23 // 
24
25
26 // --- ROOT system ---
27 #include <TClonesArray.h>
28 #include "TGeoManager.h"
29 #include "TGeoMatrix.h"
30
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliEMCALReconstructor.h"
35
36 #include "AliCodeTimer.h"
37 #include "AliCaloCalibPedestal.h"
38 #include "AliEMCALCalibData.h"
39 #include "AliESDEvent.h"
40 #include "AliESDCaloCluster.h"
41 #include "AliESDCaloCells.h"
42 #include "AliESDtrack.h"
43 #include "AliEMCALLoader.h"
44 #include "AliEMCALRawUtils.h"
45 #include "AliEMCALDigit.h"
46 #include "AliEMCALClusterizerv1.h"
47 #include "AliEMCALClusterizerv2.h"
48 #include "AliEMCALClusterizerNxN.h"
49 #include "AliEMCALRecPoint.h"
50 #include "AliEMCALPID.h"
51 #include "AliEMCALRecoUtils.h"
52 #include "AliRawReader.h"
53 #include "AliCDBEntry.h"
54 #include "AliCDBManager.h"
55 #include "AliEMCALGeometry.h"
56 #include "AliEMCAL.h"
57 #include "AliESDVZERO.h"
58 #include "AliCDBManager.h"
59 #include "AliRunLoader.h"
60 #include "AliRun.h"
61 #include "AliEMCALTriggerData.h"
62 #include "AliEMCALTriggerElectronics.h"
63 #include "AliEMCALTriggerDCSConfigDB.h"
64 #include "AliEMCALTriggerDCSConfig.h"
65 #include "AliEMCALTriggerData.h"
66 #include "AliEMCALTriggerRawDigit.h"
67 #include "AliEMCALTriggerPatch.h"
68 #include "AliEMCALTriggerTypes.h"
69
70 ClassImp(AliEMCALReconstructor) 
71   
72 const AliEMCALRecParam*     AliEMCALReconstructor::fgkRecParam        = 0;   // EMCAL rec. parameters
73 AliEMCALRawUtils*           AliEMCALReconstructor::fgRawUtils         = 0;   // EMCAL raw utilities class
74 AliEMCALClusterizer*        AliEMCALReconstructor::fgClusterizer      = 0;   // EMCAL clusterizer class
75 TClonesArray*               AliEMCALReconstructor::fgDigitsArr        = 0;   // list of digits, to be used multiple times
76 TObjArray*                  AliEMCALReconstructor::fgClustersArr      = 0;   // list of clusters, to be used multiple times
77 TClonesArray*               AliEMCALReconstructor::fgTriggerDigits    = 0;   // list of trigger digits, to be used multiple times
78 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
79 //____________________________________________________________________________
80 AliEMCALReconstructor::AliEMCALReconstructor() 
81   : fGeom(0),fCalibData(0),fPedestalData(0),fTriggerData(0x0), fMatches(0x0)
82 {
83   // ctor
84
85   // AliDebug(2, "Mark.");  
86
87   fgRawUtils = new AliEMCALRawUtils;
88   
89   //To make sure we match with the geometry in a simulation file,
90   //let's try to get it first.  If not, take the default geometry
91   AliRunLoader *rl = AliRunLoader::Instance();
92   if (rl->GetAliRun()){
93     AliEMCAL * emcal = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"));
94     if(emcal) fGeom = emcal->GetGeometry();
95   }
96   
97   if(!fGeom) {
98     AliInfo(Form("Using default geometry in reconstruction"));
99     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
100   }
101   
102   //Get calibration parameters  
103   if(!fCalibData)
104     {
105       AliCDBEntry *entry = (AliCDBEntry*) 
106         AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
107       if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
108     }
109   
110   if(!fCalibData)
111     AliFatal("Calibration parameters not found in CDB!");
112   
113   //Get calibration parameters  
114   if(!fPedestalData)
115     {
116       AliCDBEntry *entry = (AliCDBEntry*) 
117         AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
118       if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
119     }
120   
121   if(!fPedestalData)
122     AliFatal("Dead map not found in CDB!");
123   
124   if(!fGeom) AliFatal(Form("Could not get geometry!"));
125   
126   AliEMCALTriggerDCSConfigDB* dcsConfigDB = AliEMCALTriggerDCSConfigDB::Instance();
127   
128   const AliEMCALTriggerDCSConfig* dcsConfig = dcsConfigDB->GetTriggerDCSConfig();
129   
130   if (!dcsConfig) AliFatal("No Trigger DCS Configuration from OCDB!");
131   fgTriggerProcessor = new AliEMCALTriggerElectronics( dcsConfig );
132   
133   fTriggerData = new AliEMCALTriggerData();
134   
135   //Init temporary list of digits
136   fgDigitsArr     = new TClonesArray("AliEMCALDigit",1000);
137   fgClustersArr   = new TObjArray(1000);
138   fgTriggerDigits = new TClonesArray("AliEMCALTriggerRawDigit", 32 * 96);       
139         
140   //Track matching
141   fMatches = new TList();
142   fMatches->SetOwner(kTRUE);
143
144
145 //____________________________________________________________________________
146 AliEMCALReconstructor::~AliEMCALReconstructor()
147 {
148   // dtor
149
150   //AliDebug(2, "Mark.");
151
152   if(fGeom)              delete fGeom;
153   
154   //No need to delete, recovered from OCDB
155   //if(fCalibData)         delete fCalibData;
156   //if(fPedestalData)      delete fPedestalData;
157   
158   if(fgDigitsArr){
159     fgDigitsArr->Clear("C");
160     delete fgDigitsArr; 
161   }
162   
163   if(fgClustersArr){
164     fgClustersArr->Clear();
165     delete fgClustersArr; 
166   }
167   
168   if(fgTriggerDigits){
169     fgTriggerDigits->Clear();
170     delete fgTriggerDigits; 
171   }
172   
173   if(fgRawUtils)         delete fgRawUtils;
174   if(fgClusterizer)      delete fgClusterizer;
175   if(fgTriggerProcessor) delete fgTriggerProcessor;
176   
177   if(fMatches) { fMatches->Delete(); delete fMatches; fMatches = 0;}
178   
179   AliCodeTimer::Instance()->Print();
180
181
182 //____________________________________________________________________________                                  
183 void AliEMCALReconstructor::InitClusterizer() const
184 {
185   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.                          
186   Int_t clusterizerType = -1;
187   Int_t eventType = -1;
188   if(GetRecParam()) {
189     clusterizerType = GetRecParam()->GetClusterizerFlag();
190     eventType       = GetRecParam()->GetEventSpecie();
191   }
192   else{
193     AliCDBEntry *entry = (AliCDBEntry*)
194       AliCDBManager::Instance()->Get("EMCAL/Calib/RecoParam");
195     //Get The reco param for the default event specie                                                           
196     if (entry) {
197       AliEMCALRecParam *recParam  = (AliEMCALRecParam*)((TObjArray *) entry->GetObject())->At(0);
198       if(recParam) clusterizerType = recParam->GetClusterizerFlag(); 
199     }
200   }
201   
202   //Check if clusterizer previously set corresponds to what is needed for this event type                       
203   if(fgClusterizer){
204     if(eventType!=AliRecoParam::kCalib){
205       //printf("ReCreate clusterizer? Clusterizer set <%d>, Clusterizer in use <%s>\n",
206       //     clusterizerType, fgClusterizer->Version());
207       
208       if     (clusterizerType == AliEMCALRecParam::kClusterizerv1  && !strcmp(fgClusterizer->Version(),"clu-v1"))  return;
209       
210       else if(clusterizerType == AliEMCALRecParam::kClusterizerNxN && !strcmp(fgClusterizer->Version(),"clu-NxN")) return;
211       
212       else if(clusterizerType == AliEMCALRecParam::kClusterizerv2  && !strcmp(fgClusterizer->Version(),"clu-v2"))  return;
213       
214       //Need to create new clusterizer, the one set previously is not the correct one     
215       delete fgClusterizer;
216     }
217     else return;
218   }
219   
220   if      (clusterizerType  == AliEMCALRecParam::kClusterizerv1)
221     {
222       fgClusterizer = new AliEMCALClusterizerv1 (fGeom, fCalibData,fPedestalData);
223     }
224   else if (clusterizerType  == AliEMCALRecParam::kClusterizerNxN)
225     {
226       fgClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData,fPedestalData);
227     }
228   else if (clusterizerType  == AliEMCALRecParam::kClusterizerv2)
229   {
230     fgClusterizer = new AliEMCALClusterizerv2   (fGeom, fCalibData,fPedestalData);
231   }
232   else 
233   {
234     AliFatal(Form("Unknown clusterizer %d ", clusterizerType));
235   }
236 }
237
238 //____________________________________________________________________________
239 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
240 {
241   // method called by AliReconstruction; 
242   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
243   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
244   // the global tracking.
245   // Works on the current event.
246   
247   AliCodeTimerAuto("",0)
248     
249   //Get input digits and put them in fgDigitsArr, clear the list before 
250   ReadDigitsArrayFromTree(digitsTree);
251   
252   InitClusterizer();
253   
254   fgClusterizer->InitParameters();
255   fgClusterizer->SetOutput(clustersTree);
256   
257   //Skip clusterization of LED events
258   if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
259     
260     if(fgDigitsArr && fgDigitsArr->GetEntries()) {
261       
262       fgClusterizer->SetInput(digitsTree);
263       
264       //fgClusterizer->Digits2Clusters("deb all") ; //For debugging
265       fgClusterizer->Digits2Clusters("");
266       
267       fgClusterizer->Clear();
268       
269     }//digits array exists and has somethind
270   }//not a LED event
271   
272   clustersTree->Fill(); 
273
274   // Deleting the recpoints at the end of the reconstruction call
275   fgClusterizer->DeleteRecPoints();
276 }
277
278 //____________________________________________________________________________
279 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
280   
281 {
282   // Conversion from raw data to
283   // EMCAL digits.
284   // Works on a single-event basis
285   
286   rawReader->Reset() ; 
287   
288   fTriggerData->SetMode(1);     
289   
290   if(fgDigitsArr) fgDigitsArr->Clear("C");
291   
292   TClonesArray *digitsTrg = new TClonesArray("AliEMCALTriggerRawDigit", 32 * 96);
293   
294   Int_t bufsize = 32000;
295   digitsTree->Branch("EMCAL", &fgDigitsArr, bufsize);
296   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
297   
298   //Skip calibration events do the rest
299   Bool_t doFit = kTRUE;
300   if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
301   if (doFit){
302     //must be done here because, in constructor, option is not yet known
303     fgRawUtils->SetOption(GetOption());
304     
305     //  fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
306     
307     //   fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
308     //    fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
309     fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
310     fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
311     fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
312     if (!fgRawUtils->GetFittingAlgorithm()) fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
313     fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
314     
315     //fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
316     //fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
317     
318     //  fgRawUtils->SetTimeMin(-99999 );
319     //  fgRawUtils->SetTimeMax( 99999 );
320     
321     fgRawUtils->Raw2Digits(rawReader,fgDigitsArr,fPedestalData,digitsTrg,fTriggerData);
322     
323   }//skip calibration event
324   else{
325     AliDebug(1," Calibration Event, skip!");
326   }
327   
328   digitsTree->Fill();
329   digitsTrg->Delete();
330   delete digitsTrg;
331   
332 }
333
334
335 //____________________________________________________________________________
336 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
337                                     AliESDEvent* esd) const
338 {
339   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
340   // and V0 
341   // Works on the current event
342   // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
343   //return;
344   
345   //########################################
346   // Trigger
347   //########################################
348   
349   static int saveOnce = 0;
350   
351   Int_t v0M[2] = {0, 0};
352   
353   AliESDVZERO* esdV0 = esd->GetVZEROData();
354   
355   if (esdV0) 
356     {
357                 v0M[0] = esdV0->GetTriggerChargeC();
358                 v0M[1] = esdV0->GetTriggerChargeA();
359     }
360   else
361     {
362       AliWarning("No V0 ESD! Run trigger processor w/ null V0 charges");
363     }
364   
365   if (fgTriggerDigits && fgTriggerDigits->GetEntriesFast()) fgTriggerDigits->Delete();
366   
367   TBranch *branchtrg = digitsTree->GetBranch("EMTRG");
368   
369   if (!branchtrg) 
370     { 
371       AliError("Can't get the branch with the EMCAL trigger digits!");
372       return;
373     }
374   
375   branchtrg->SetAddress(&fgTriggerDigits);
376   branchtrg->GetEntry(0);
377   
378   // Note: fgTriggerProcessor reset done at the end of this method
379   fgTriggerProcessor->Digits2Trigger(fgTriggerDigits, v0M, fTriggerData);
380   
381   // Fill ESD
382   AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("EMCAL");
383   
384   if (trgESD)
385     {
386       trgESD->Allocate(fgTriggerDigits->GetEntriesFast());
387       
388       for (Int_t i = 0; i < fgTriggerDigits->GetEntriesFast(); i++)
389         {         
390           AliEMCALTriggerRawDigit* rdig = (AliEMCALTriggerRawDigit*)fgTriggerDigits->At(i);
391           if (AliDebugLevel() > 999) rdig->Print("");
392                 
393           Int_t px, py;
394           if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
395             {
396               Int_t a = -1, t = -1, times[10]; 
397               
398               rdig->GetMaximum(a, t);
399               rdig->GetL0Times(times);
400                         
401               trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum(), rdig->GetTriggerBits());
402             }
403         }
404       
405                 for (int i = 0; i < 2; i++) {
406                         trgESD->SetL1Threshold(2 * i    , fTriggerData->GetL1JetThreshold(  i));
407                         trgESD->SetL1Threshold(2 * i + 1, fTriggerData->GetL1GammaThreshold(i));
408                 }
409       
410       Int_t v0[2];
411       fTriggerData->GetL1V0(v0);
412       
413       trgESD->SetL1V0(v0);      
414       trgESD->SetL1FrameMask(fTriggerData->GetL1FrameMask());            
415       
416       if (!saveOnce && fTriggerData->GetL1DataDecoded()) 
417         {
418           int type[15] = {0};
419           fTriggerData->GetL1TriggerType(type);
420           
421           esd->SetCaloTriggerType(type);
422           
423           saveOnce = 1;
424         }
425     }
426   
427   // Resetting
428   fTriggerData->Reset();
429   
430   //########################################
431   //##############Fill CaloCells###############
432   //########################################
433   
434   //Get input digits and put them in fgDigitsArr, clear the list before 
435   ReadDigitsArrayFromTree(digitsTree);
436   
437   Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
438   AliDebug(1,Form("%d digits",nDigits));
439   
440   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
441   emcCells.CreateContainer(nDigits);
442   emcCells.SetType(AliVCaloCells::kEMCALCell);
443   
444   Float_t energy = 0;
445   Float_t time   = 0;
446   for (Int_t idig = 0 ; idig < nDigits ; idig++) 
447   {
448     const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
449     time   = dig->GetTime();      // Time already calibrated in clusterizer
450     energy = dig->GetAmplitude(); // energy calibrated in clusterizer
451     
452     if(energy > 0 )
453     {
454       fgClusterizer->Calibrate(energy,time,dig->GetId()); //Digits already calibrated in clusterizers
455       
456       if(energy > 0) //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
457       { 
458         // Only for MC
459         // Get the label of the primary particle that generated the cell
460         // Assign the particle that deposited more energy
461         Int_t   nprimaries = dig->GetNprimary() ;
462         Int_t   digLabel   =-1 ;
463         Float_t edep       =-1.;
464         if ( nprimaries > 0 )
465         {
466           Int_t jndex ;
467           for ( jndex = 0 ; jndex < nprimaries ; jndex++ ) { // all primaries in digit
468  
469             if(edep < dig->GetDEParent(jndex+1))
470             {
471               digLabel = dig->GetIparent (jndex+1);
472               edep     = dig->GetDEParent(jndex+1);
473             }
474                    
475           } // all primaries in digit      
476         } // select primary label
477         
478         emcCells.SetCell(idignew,dig->GetId(),energy, time,digLabel);   
479         idignew++;
480       }
481     }
482   }
483   
484   emcCells.SetNumberOfCells(idignew);
485   emcCells.Sort();
486   
487   //------------------------------------------------------------
488   //-----------------CLUSTERS-----------------------------
489   //------------------------------------------------------------
490   clustersTree->SetBranchStatus("*",0); //disable all branches
491   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
492   if(fgClustersArr) fgClustersArr->Clear();
493   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
494   branch->SetAddress(&fgClustersArr);
495   branch->GetEntry(0);
496   //clustersTree->GetEvent(0);
497   
498   Int_t nClusters = fgClustersArr->GetEntries(),  nClustersNew=0;
499   AliDebug(1,Form("%d clusters",nClusters));
500
501   
502   //########################################
503   //##############Fill CaloClusters#############
504   //########################################
505   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
506     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
507     if(!clust) continue;
508     //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
509     // clust->Print(); //For debugging
510     // Get information from EMCAL reconstruction points
511     Float_t xyz[3];
512     TVector3 gpos;
513     clust->GetGlobalPosition(gpos);
514     for (Int_t ixyz=0; ixyz<3; ixyz++)
515       xyz[ixyz] = gpos[ixyz];
516     Float_t elipAxis[2];
517     clust->GetElipsAxis(elipAxis);
518     //Create digits lists
519     Int_t cellMult = clust->GetMultiplicity();
520     //TArrayS digiList(digitMult);
521     Float_t *amplFloat = clust->GetEnergiesList();
522     Int_t   *digitInts = clust->GetAbsId();
523     TArrayS absIdList(cellMult);
524     TArrayD fracList(cellMult);
525     
526     Int_t newCellMult = 0;
527     for (Int_t iCell=0; iCell<cellMult; iCell++) {
528       if (amplFloat[iCell] > 0) {
529         absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
530         //Calculate Fraction
531         if(emcCells.GetCellAmplitude(digitInts[iCell])>0 && GetRecParam()->GetUnfold()){
532           fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell]));//get cell calibration value 
533           
534         }
535         else{
536           fracList[newCellMult] = 0; 
537         }
538         newCellMult++;
539       }
540     }
541     
542     absIdList.Set(newCellMult);
543     fracList.Set(newCellMult);
544     
545     if(newCellMult > 0) { // accept cluster if it has some digit
546       nClustersNew++;
547       //Primaries
548       Int_t  parentMult  = 0;
549       Int_t *parentList =  clust->GetParents(parentMult);
550       // fills the ESDCaloCluster
551       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
552       ec->SetType(AliVCluster::kEMCALClusterv1);
553       ec->SetPosition(xyz);
554       ec->SetE(clust->GetEnergy());
555       
556       //Distance to the nearest bad crystal
557       ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
558       
559       ec->SetNCells(newCellMult);
560       //Change type of list from short to ushort
561       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
562       Double_t *newFracList   = new Double_t[newCellMult];
563       for(Int_t i = 0; i < newCellMult ; i++) {
564         newAbsIdList[i]=absIdList[i];
565         newFracList[i] =fracList[i];
566       }
567       ec->SetCellsAbsId(newAbsIdList);
568       ec->SetCellsAmplitudeFraction(newFracList);
569       ec->SetDispersion(clust->GetDispersion());
570       ec->SetChi2(-1); //not yet implemented
571       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
572       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
573       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
574       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
575   
576       
577       TArrayI arrayParents(parentMult,parentList);
578       ec->AddLabels(arrayParents);
579       //
580       //Track matching
581       //
582       fMatches->Clear();
583       Int_t nTracks = esd->GetNumberOfTracks();
584       for (Int_t itrack = 0; itrack < nTracks; itrack++)
585         {
586           AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
587           if(track->GetEMCALcluster()==iClust)
588             {
589               Float_t dEta=-999, dPhi=-999;
590               Bool_t isMatch =  CalculateResidual(track, ec, dEta, dPhi);
591               if(!isMatch) 
592                 {
593                   // AliDebug(10, "Not good");
594                   continue;
595                 }
596               AliEMCALMatch *match = new AliEMCALMatch();
597               match->SetIndexT(itrack);
598               match->SetDistance(TMath::Sqrt(dEta*dEta+dPhi*dPhi));
599               match->SetdEta(dEta);
600               match->SetdPhi(dPhi);
601               fMatches->Add(match);
602             }
603         } 
604       fMatches->Sort(kSortAscending); //Sort matched tracks from closest to furthest
605       Int_t nMatch = fMatches->GetEntries();
606       TArrayI arrayTrackMatched(nMatch);
607       for(Int_t imatch=0; imatch<nMatch; imatch++)
608         {
609           AliEMCALMatch *match = (AliEMCALMatch*)fMatches->At(imatch);
610           arrayTrackMatched[imatch] = match->GetIndexT();
611           if(imatch==0)
612             {
613               ec->SetTrackDistance(match->GetdPhi(), match->GetdEta());
614             }
615         }
616       ec->AddTracksMatched(arrayTrackMatched);
617     
618       //add the cluster to the esd object
619       esd->AddCaloCluster(ec);
620
621       delete ec;
622       delete [] newAbsIdList ;
623       delete [] newFracList ;
624     }
625   } // cycle on clusters
626
627   //
628   //Reset the index of matched cluster for tracks
629   //to the one in CaloCluster array
630   Int_t ncls = esd->GetNumberOfCaloClusters();
631   for(Int_t icl=0; icl<ncls; icl++)
632     {
633       AliESDCaloCluster *cluster = esd->GetCaloCluster(icl);
634       if(!cluster || !cluster->IsEMCAL()) continue;
635       TArrayI *trackIndex = cluster->GetTracksMatched();
636       for(Int_t itr=0; itr<trackIndex->GetSize(); itr++)
637         {
638           AliESDtrack *track = esd->GetTrack(trackIndex->At(itr));
639           track->SetEMCALcluster(cluster->GetID());
640         }
641     }
642   
643   
644   //Fill ESDCaloCluster with PID weights
645   AliEMCALPID *pid = new AliEMCALPID;
646   //pid->SetPrintInfo(kTRUE);
647   pid->SetReconstructor(kTRUE);
648   pid->RunPID(esd);
649   delete pid;
650   
651   //Store EMCAL misalignment matrixes
652   FillMisalMatrixes(esd) ;
653   
654 }
655
656 //==================================================================================
657 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
658   //Store EMCAL matrixes in ESD Header
659   
660   //Check, if matrixes was already stored
661   for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
662     if(esd->GetEMCALMatrix(sm)!=0)
663       return ;
664   }
665   
666   //Create and store matrixes
667   if(!gGeoManager){
668     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
669     return ;
670   }
671   //Note, that owner of copied marixes will be header
672   const Int_t bufsize = 255;
673   char path[bufsize] ;
674   TGeoHMatrix * m = 0x0;
675   for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
676     snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
677     if(sm >= 10 && !((fGeom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
678     if(sm >= 10 &&  ((fGeom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM3rd_%d",sm-10+1) ;
679     
680     if (gGeoManager->CheckPath(path)){
681       gGeoManager->cd(path);
682       m = gGeoManager->GetCurrentMatrix() ;
683       //                        printf("================================================= \n");
684       //                        printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
685       //                        m->Print("");
686       esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
687       //                        printf("================================================= \n");
688     }
689     else{
690       esd->SetEMCALMatrix(NULL,sm) ;
691     }
692   }
693 }
694
695 //__________________________________________________________________________
696 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
697 {
698   // Read the digits from the input tree
699   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);    
700   
701   // Clear previous digits in the list
702   if(fgDigitsArr){ 
703     fgDigitsArr->Clear("C");
704   }
705   else{
706     // It should not happen, but just in case ...
707     fgDigitsArr = new TClonesArray("AliEMCALDigit",100); 
708   }
709   
710   // Read the digits from the input tree
711   TBranch *branch = digitsTree->GetBranch("EMCAL");
712   if (!branch) { 
713     AliError("can't get the branch with the EMCAL digits !");
714     return;
715   }  
716   
717   branch->SetAddress(&fgDigitsArr);
718   branch->GetEntry(0);
719 }
720
721 //==================================================================================
722 Bool_t AliEMCALReconstructor::CalculateResidual(AliESDtrack *track, AliESDCaloCluster *cluster, Float_t &dEta, Float_t &dPhi)const
723 {
724   //
725   // calculate the residual between track and cluster
726   //
727
728   // If the esdFriend is available, use the TPCOuter point as the starting point of extrapolation
729   // Otherwise use the TPCInner point
730
731   dEta = -999, dPhi = -999;
732
733   AliExternalTrackParam *trkParam = 0;
734   const AliESDfriendTrack*  friendTrack = track->GetFriendTrack();
735   if(friendTrack && friendTrack->GetTPCOut())
736     trkParam = const_cast<AliExternalTrackParam*>(friendTrack->GetTPCOut());
737   else
738     trkParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
739   if(!trkParam) return kFALSE;
740
741   AliExternalTrackParam trkParamTmp (*trkParam);
742   if(!AliEMCALRecoUtils::ExtrapolateTrackToCluster(&trkParamTmp, cluster, track->GetMass(kTRUE), GetRecParam()->GetExtrapolateStep(), dEta, dPhi)) return kFALSE;
743
744   return kTRUE;
745 }
746
747 //
748 //==================================================================================
749 //
750 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch() 
751   : TObject(),  
752     fIndexT(-1), 
753     fDistance(-999.),
754     fdEta(-999.),
755     fdPhi(-999.)
756 {
757   //default constructor
758
759 }
760
761 //
762 //==================================================================================
763 //
764 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch(const AliEMCALMatch& copy)
765   : TObject(),
766     fIndexT(copy.fIndexT),
767     fDistance(copy.fDistance),
768     fdEta(copy.fdEta),
769     fdPhi(copy.fdPhi)
770 {
771   //copy ctor
772 }
773 //_____________________________________________________________________
774 AliEMCALReconstructor::AliEMCALMatch& AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch::operator = (const AliEMCALMatch &source)
775 { // assignment operator; use copy ctor
776   if (&source == this) return *this;
777
778   new (this) AliEMCALMatch(source);
779   return *this;
780 }
781 //
782 //==================================================================================
783 //
784 Int_t AliEMCALReconstructor::AliEMCALMatch::Compare(const TObject *obj) const 
785 {
786   //
787   // Compare wrt the residual
788   //
789         
790   AliEMCALReconstructor::AliEMCALMatch *that = (AliEMCALReconstructor::AliEMCALMatch*)obj;
791         
792   Double_t thisDist = fDistance;//fDistance;
793   Double_t thatDist = that->fDistance;//that->GetDistance();
794         
795   if (thisDist > thatDist) return 1;
796   else if (thisDist < thatDist) return -1;
797   return 0;
798 }