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