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