]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Bug #92237 fixed as Barth suggested
[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",1000);   
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->Clear();
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           
392           Int_t px, py;
393           if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
394             {
395               Int_t a = -1, t = -1, times[10]; 
396               
397               rdig->GetMaximum(a, t);
398               rdig->GetL0Times(times);
399               
400               trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum(), rdig->GetTriggerBits());
401             }
402         }
403       
404       trgESD->SetL1Threshold(0, fTriggerData->GetL1GammaThreshold());
405       
406       trgESD->SetL1Threshold(1, fTriggerData->GetL1JetThreshold()  );
407       
408       Int_t v0[2];
409       fTriggerData->GetL1V0(v0);
410       
411       trgESD->SetL1V0(v0);      
412       trgESD->SetL1FrameMask(fTriggerData->GetL1FrameMask());            
413       
414       if (!saveOnce && fTriggerData->GetL1DataDecoded()) 
415         {
416           int type[8] = {0};
417           fTriggerData->GetL1TriggerType(type);
418           
419           esd->SetCaloTriggerType(type);
420           
421           saveOnce = 1;
422         }
423     }
424   
425   // Resetting
426   fTriggerData->Reset();
427   
428   //########################################
429   //##############Fill CaloCells###############
430   //########################################
431   
432   //Get input digits and put them in fgDigitsArr, clear the list before 
433   ReadDigitsArrayFromTree(digitsTree);
434   
435   Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
436   AliDebug(1,Form("%d digits",nDigits));
437   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
438   emcCells.CreateContainer(nDigits);
439   emcCells.SetType(AliVCaloCells::kEMCALCell);
440   Float_t energy = 0;
441   Float_t time   = 0;
442   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
443     const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
444     time   = dig->GetTime();      // Time already calibrated in clusterizer
445     energy = dig->GetAmplitude(); // energy calibrated in clusterizer
446     if(energy > 0 ){
447       fgClusterizer->Calibrate(energy,time,dig->GetId()); //Digits already calibrated in clusterizers
448       if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
449         emcCells.SetCell(idignew,dig->GetId(),energy, time);   
450         idignew++;
451       }
452     }
453   }
454   emcCells.SetNumberOfCells(idignew);
455   emcCells.Sort();
456   
457   //------------------------------------------------------------
458   //-----------------CLUSTERS-----------------------------
459   //------------------------------------------------------------
460   clustersTree->SetBranchStatus("*",0); //disable all branches
461   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
462   if(fgClustersArr) fgClustersArr->Clear();
463   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
464   branch->SetAddress(&fgClustersArr);
465   branch->GetEntry(0);
466   //clustersTree->GetEvent(0);
467   
468   Int_t nClusters = fgClustersArr->GetEntries(),  nClustersNew=0;
469   AliDebug(1,Form("%d clusters",nClusters));
470
471   
472   //########################################
473   //##############Fill CaloClusters#############
474   //########################################
475   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
476     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
477     if(!clust) continue;
478     //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
479     // clust->Print(); //For debugging
480     // Get information from EMCAL reconstruction points
481     Float_t xyz[3];
482     TVector3 gpos;
483     clust->GetGlobalPosition(gpos);
484     for (Int_t ixyz=0; ixyz<3; ixyz++)
485       xyz[ixyz] = gpos[ixyz];
486     Float_t elipAxis[2];
487     clust->GetElipsAxis(elipAxis);
488     //Create digits lists
489     Int_t cellMult = clust->GetMultiplicity();
490     //TArrayS digiList(digitMult);
491     Float_t *amplFloat = clust->GetEnergiesList();
492     Int_t   *digitInts = clust->GetAbsId();
493     TArrayS absIdList(cellMult);
494     TArrayD fracList(cellMult);
495     
496     Int_t newCellMult = 0;
497     for (Int_t iCell=0; iCell<cellMult; iCell++) {
498       if (amplFloat[iCell] > 0) {
499         absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
500         //Calculate Fraction
501         if(emcCells.GetCellAmplitude(digitInts[iCell])>0 && GetRecParam()->GetUnfold()){
502           fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell]));//get cell calibration value 
503           
504         }
505         else{
506           fracList[newCellMult] = 0; 
507         }
508         newCellMult++;
509       }
510     }
511     
512     absIdList.Set(newCellMult);
513     fracList.Set(newCellMult);
514     
515     if(newCellMult > 0) { // accept cluster if it has some digit
516       nClustersNew++;
517       //Primaries
518       Int_t  parentMult  = 0;
519       Int_t *parentList =  clust->GetParents(parentMult);
520       // fills the ESDCaloCluster
521       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
522       ec->SetType(AliVCluster::kEMCALClusterv1);
523       ec->SetPosition(xyz);
524       ec->SetE(clust->GetEnergy());
525       
526       //Distance to the nearest bad crystal
527       ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
528       
529       ec->SetNCells(newCellMult);
530       //Change type of list from short to ushort
531       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
532       Double_t *newFracList   = new Double_t[newCellMult];
533       for(Int_t i = 0; i < newCellMult ; i++) {
534         newAbsIdList[i]=absIdList[i];
535         newFracList[i] =fracList[i];
536       }
537       ec->SetCellsAbsId(newAbsIdList);
538       ec->SetCellsAmplitudeFraction(newFracList);
539       ec->SetDispersion(clust->GetDispersion());
540       ec->SetChi2(-1); //not yet implemented
541       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
542       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
543       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
544       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
545   
546       
547       TArrayI arrayParents(parentMult,parentList);
548       ec->AddLabels(arrayParents);
549       //
550       //Track matching
551       //
552       fMatches->Clear();
553       Int_t nTracks = esd->GetNumberOfTracks();
554       for (Int_t itrack = 0; itrack < nTracks; itrack++)
555         {
556           AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
557           if(track->GetEMCALcluster()==iClust)
558             {
559               Float_t dEta=-999, dPhi=-999;
560               Bool_t isMatch =  CalculateResidual(track, ec, dEta, dPhi);
561               if(!isMatch) 
562                 {
563                   // AliDebug(10, "Not good");
564                   continue;
565                 }
566               AliEMCALMatch *match = new AliEMCALMatch();
567               match->SetIndexT(itrack);
568               match->SetDistance(TMath::Sqrt(dEta*dEta+dPhi*dPhi));
569               match->SetdEta(dEta);
570               match->SetdPhi(dPhi);
571               fMatches->Add(match);
572             }
573         } 
574       fMatches->Sort(kSortAscending); //Sort matched tracks from closest to furthest
575       Int_t nMatch = fMatches->GetEntries();
576       TArrayI arrayTrackMatched(nMatch);
577       for(Int_t imatch=0; imatch<nMatch; imatch++)
578         {
579           AliEMCALMatch *match = (AliEMCALMatch*)fMatches->At(imatch);
580           arrayTrackMatched[imatch] = match->GetIndexT();
581           if(imatch==0)
582             {
583               ec->SetTrackDistance(match->GetdPhi(), match->GetdEta());
584             }
585         }
586       ec->AddTracksMatched(arrayTrackMatched);
587     
588       //add the cluster to the esd object
589       esd->AddCaloCluster(ec);
590
591       delete ec;
592       delete [] newAbsIdList ;
593       delete [] newFracList ;
594     }
595   } // cycle on clusters
596
597   //
598   //Reset the index of matched cluster for tracks
599   //to the one in CaloCluster array
600   Int_t ncls = esd->GetNumberOfCaloClusters();
601   for(Int_t icl=0; icl<ncls; icl++)
602     {
603       AliESDCaloCluster *cluster = esd->GetCaloCluster(icl);
604       if(!cluster || !cluster->IsEMCAL()) continue;
605       TArrayI *trackIndex = cluster->GetTracksMatched();
606       for(Int_t itr=0; itr<trackIndex->GetSize(); itr++)
607         {
608           AliESDtrack *track = esd->GetTrack(trackIndex->At(itr));
609           track->SetEMCALcluster(cluster->GetID());
610         }
611     }
612   
613   
614   //Fill ESDCaloCluster with PID weights
615   AliEMCALPID *pid = new AliEMCALPID;
616   //pid->SetPrintInfo(kTRUE);
617   pid->SetReconstructor(kTRUE);
618   pid->RunPID(esd);
619   delete pid;
620   
621   //Store EMCAL misalignment matrixes
622   FillMisalMatrixes(esd) ;
623   
624 }
625
626 //==================================================================================
627 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
628   //Store EMCAL matrixes in ESD Header
629   
630   //Check, if matrixes was already stored
631   for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
632     if(esd->GetEMCALMatrix(sm)!=0)
633       return ;
634   }
635   
636   //Create and store matrixes
637   if(!gGeoManager){
638     AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
639     return ;
640   }
641   //Note, that owner of copied marixes will be header
642   const Int_t bufsize = 255;
643   char path[bufsize] ;
644   TGeoHMatrix * m = 0x0;
645   for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
646     snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
647     if(sm >= 10 && !((fGeom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
648     if(sm >= 10 &&  ((fGeom->GetEMCGeometry()->GetGeoName()).Contains("12SMV1"))) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM3rd_%d",sm-10+1) ;
649     
650     if (gGeoManager->CheckPath(path)){
651       gGeoManager->cd(path);
652       m = gGeoManager->GetCurrentMatrix() ;
653       //                        printf("================================================= \n");
654       //                        printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
655       //                        m->Print("");
656       esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
657       //                        printf("================================================= \n");
658     }
659     else{
660       esd->SetEMCALMatrix(NULL,sm) ;
661     }
662   }
663 }
664
665 //__________________________________________________________________________
666 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
667 {
668   // Read the digits from the input tree
669   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);    
670   
671   // Clear previous digits in the list
672   if(fgDigitsArr){ 
673     fgDigitsArr->Clear("C");
674   }
675   else{
676     // It should not happen, but just in case ...
677     fgDigitsArr = new TClonesArray("AliEMCALDigit",100); 
678   }
679   
680   // Read the digits from the input tree
681   TBranch *branch = digitsTree->GetBranch("EMCAL");
682   if (!branch) { 
683     AliError("can't get the branch with the EMCAL digits !");
684     return;
685   }  
686   
687   branch->SetAddress(&fgDigitsArr);
688   branch->GetEntry(0);
689 }
690
691 //==================================================================================
692 Bool_t AliEMCALReconstructor::CalculateResidual(AliESDtrack *track, AliESDCaloCluster *cluster, Float_t &dEta, Float_t &dPhi)const
693 {
694   //
695   // calculate the residual between track and cluster
696   //
697
698   // If the esdFriend is available, use the TPCOuter point as the starting point of extrapolation
699   // Otherwise use the TPCInner point
700
701   dEta = -999, dPhi = -999;
702
703   AliExternalTrackParam *trkParam = 0;
704   const AliESDfriendTrack*  friendTrack = track->GetFriendTrack();
705   if(friendTrack && friendTrack->GetTPCOut())
706     trkParam = const_cast<AliExternalTrackParam*>(friendTrack->GetTPCOut());
707   else
708     trkParam = const_cast<AliExternalTrackParam*>(track->GetInnerParam());
709   if(!trkParam) return kFALSE;
710
711   AliExternalTrackParam trkParamTmp (*trkParam);
712   if(!AliEMCALRecoUtils::ExtrapolateTrackToCluster(&trkParamTmp, cluster, track->GetMass(), GetRecParam()->GetExtrapolateStep(), dEta, dPhi)) return kFALSE;
713
714   return kTRUE;
715 }
716
717 //
718 //==================================================================================
719 //
720 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch() 
721   : TObject(),  
722     fIndexT(-1), 
723     fDistance(-999.),
724     fdEta(-999.),
725     fdPhi(-999.)
726 {
727   //default constructor
728
729 }
730
731 //
732 //==================================================================================
733 //
734 AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch(const AliEMCALMatch& copy)
735   : TObject(),
736     fIndexT(copy.fIndexT),
737     fDistance(copy.fDistance),
738     fdEta(copy.fdEta),
739     fdPhi(copy.fdPhi)
740 {
741   //copy ctor
742 }
743 //_____________________________________________________________________
744 AliEMCALReconstructor::AliEMCALMatch& AliEMCALReconstructor::AliEMCALMatch::AliEMCALMatch::operator = (const AliEMCALMatch &source)
745 { // assignment operator; use copy ctor
746   if (&source == this) return *this;
747
748   new (this) AliEMCALMatch(source);
749   return *this;
750 }
751 //
752 //==================================================================================
753 //
754 Int_t AliEMCALReconstructor::AliEMCALMatch::Compare(const TObject *obj) const 
755 {
756   //
757   // Compare wrt the residual
758   //
759         
760   AliEMCALReconstructor::AliEMCALMatch *that = (AliEMCALReconstructor::AliEMCALMatch*)obj;
761         
762   Double_t thisDist = fDistance;//fDistance;
763   Double_t thatDist = that->fDistance;//that->GetDistance();
764         
765   if (thisDist > thatDist) return 1;
766   else if (thisDist < thatDist) return -1;
767   return 0;
768 }