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