]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
remove from digits array those channels with bad timing or which only produced Low...
[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 //-- Aleksei Pavlinov : added staf for EMCAL jet trigger 9Apr 25, 2008)
25 //                    : fgDigitsArr should read just once at event
26
27 // --- ROOT system ---
28 #include <TList.h>
29 #include <TClonesArray.h>
30 #include <TH2.h>
31 #include "TGeoManager.h"
32 #include "TGeoMatrix.h"
33
34 // --- Standard library ---
35
36 // --- AliRoot header files ---
37 #include "AliEMCALReconstructor.h"
38
39 #include "AliCodeTimer.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 "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 "AliVZEROLoader.h"
63
64 ClassImp(AliEMCALReconstructor) 
65
66 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0;  // EMCAL rec. parameters
67 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0;   // EMCAL raw utilities class
68 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0;   // EMCAL clusterizer class
69 TClonesArray*     AliEMCALReconstructor::fgDigitsArr = 0;  // shoud read just once at event
70 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
71 //____________________________________________________________________________
72 AliEMCALReconstructor::AliEMCALReconstructor() 
73   : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0) 
74 {
75   // ctor
76
77   fgRawUtils = new AliEMCALRawUtils;
78
79   //To make sure we match with the geometry in a simulation file,
80   //let's try to get it first.  If not, take the default geometry
81   AliRunLoader *rl = AliRunLoader::Instance();
82   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
83     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
84   } else {
85     AliInfo(Form("Using default geometry in reconstruction"));
86     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
87   }
88
89   //Get calibration parameters  
90   if(!fCalibData)
91     {
92                 AliCDBEntry *entry = (AliCDBEntry*) 
93                 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
94                 if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
95     }
96         
97   if(!fCalibData)
98                 AliFatal("Calibration parameters not found in CDB!");
99         
100   //Get calibration parameters  
101   if(!fPedestalData)
102     {
103                 AliCDBEntry *entry = (AliCDBEntry*) 
104                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
105                 if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
106     }
107         
108         if(!fPedestalData)
109                 AliFatal("Dead map not found in CDB!");
110         
111         
112   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
113   fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData); 
114         
115   if(!fGeom) AliFatal(Form("Could not get geometry!"));
116
117   fgTriggerProcessor = new AliEMCALTriggerElectronics();
118
119
120 //____________________________________________________________________________
121 AliEMCALReconstructor::~AliEMCALReconstructor()
122 {
123   // dtor
124   delete fGeom;
125   delete fgRawUtils;
126   delete fgClusterizer;
127   delete fgTriggerProcessor;
128
129   AliCodeTimer::Instance()->Print();
130
131
132 // //____________________________________________________________________________
133 // void AliEMCALReconstructor::Init()
134 // {
135 //   // Trigger hists - Oct 24, 2007
136 //    fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
137 // }
138
139 //____________________________________________________________________________
140 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
141 {
142   // method called by AliReconstruction; 
143   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
144   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
145   // the global tracking.
146   // Works on the current event.
147
148   AliCodeTimerAuto("",0)
149
150   ReadDigitsArrayFromTree(digitsTree);
151   fgClusterizer->InitParameters();
152   fgClusterizer->SetOutput(clustersTree);
153
154   AliEMCALTriggerData* trgData = new AliEMCALTriggerData();
155         
156   Int_t bufferSize = 32000;
157
158   if (TBranch* triggerBranch = clustersTree->GetBranch("EMTRG"))
159           triggerBranch->SetAddress(&trgData);
160   else
161           clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
162
163 //  AliVZEROLoader* vzeroLoader = dynamic_cast<AliVZEROLoader*>(AliRunLoader::Instance()->GetDetectorLoader("VZERO"));
164 //  
165 //  TTree* treeV0 = 0x0;
166 //      
167 //  if (vzeroLoader) 
168 //  {
169 //        vzeroLoader->LoadDigits("READ");
170 //      treeV0 = vzeroLoader->TreeD();
171 //  }
172 //
173   TClonesArray *trgDigits = new TClonesArray("AliEMCALRawDigit",1000);
174   TBranch *branchdig = digitsTree->GetBranch("EMTRG");
175   if (!branchdig) 
176   { 
177           AliError("Can't get the branch with the EMCAL trigger digits !");
178           return;
179   }
180
181   branchdig->SetAddress(&trgDigits);
182   branchdig->GetEntry(0);
183
184   fgTriggerProcessor->Digits2Trigger(trgDigits, NULL, trgData);
185         
186   trgDigits->Delete();
187   delete        trgDigits;
188
189   if(fgDigitsArr && fgDigitsArr->GetEntries()) {
190
191     fgClusterizer->SetInput(digitsTree);
192     
193     if(Debug())
194       fgClusterizer->Digits2Clusters("deb all") ;
195     else
196       fgClusterizer->Digits2Clusters("");
197     
198     fgClusterizer->Clear();
199
200   }
201
202 //  if (vzeroLoader) 
203 //  {
204 //        vzeroLoader->UnloadDigits();
205 //  }
206 //
207   clustersTree->Fill(); 
208
209   delete trgData;
210 }
211
212 //____________________________________________________________________________
213 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
214
215 {
216   // Conversion from raw data to
217   // EMCAL digits.
218   // Works on a single-event basis
219
220   rawReader->Reset() ; 
221
222   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
223   TClonesArray *digitsTrg = new TClonesArray("AliEMCALRawDigit", 200);
224
225   Int_t bufsize = 32000;
226   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
227   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
228
229   //must be done here because, in constructor, option is not yet known
230   fgRawUtils->SetOption(GetOption());
231
232   fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
233   fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
234   fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
235   fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
236   fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
237   fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
238   fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
239   fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
240   fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
241   fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
242         
243   fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
244
245   digitsTree->Fill();
246   digitsArr->Delete();
247   digitsTrg->Delete();
248   delete digitsArr;
249   delete digitsTrg;
250
251 }
252
253
254 //____________________________________________________________________________
255 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
256                                     AliESDEvent* esd) const
257 {
258   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
259   // and V0 
260   // Works on the current event
261   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
262   //return;
263
264   //######################################################
265   //#########Calculate trigger and set trigger info###########
266   //######################################################
267   // Obsolete, to be changed with new trigger emulator when consensus is achieved about what is stored in ESDs.
268   AliEMCALTrigger tr;
269   //   tr.SetPatchSize(1);  // create 4x4 patches
270   tr.SetSimulation(kFALSE); // Reconstruction mode
271   tr.SetDigitsList(fgDigitsArr);
272   // Get VZERO total multiplicity for jet trigger simulation 
273   // The simulation of jey trigger will be incorrect if no VZERO data 
274   // at ESD
275   AliESDVZERO* vZero = esd->GetVZEROData();
276   if(vZero) {
277     tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
278   }
279   //
280   tr.Trigger();
281
282   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
283   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
284   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
285   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
286
287   Int_t iSM2x2      = tr.Get2x2SuperModule();
288   Int_t iSMnxn      = tr.GetnxnSuperModule();
289   Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
290   Int_t iModulePhinxn = tr.GetnxnModulePhi();
291   Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
292   Int_t iModuleEtanxn = tr.GetnxnModuleEta();
293
294   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
295   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",  maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
296
297   TVector3    pos2x2(-1,-1,-1);
298   TVector3    posnxn(-1,-1,-1);
299
300   Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
301   Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
302   fGeom->GetGlobal(iAbsId2x2, pos2x2);
303   fGeom->GetGlobal(iAbsIdnxn, posnxn);
304   //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
305   
306   TArrayF triggerPosition(6);
307   triggerPosition[0] = pos2x2(0) ;   
308   triggerPosition[1] = pos2x2(1) ;   
309   triggerPosition[2] = pos2x2(2) ;  
310   triggerPosition[3] = posnxn(0) ;   
311   triggerPosition[4] = posnxn(1) ;   
312   triggerPosition[5] = posnxn(2) ;
313   //printf(" triggerPosition ");
314   //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
315
316   TArrayF triggerAmplitudes(4);
317   triggerAmplitudes[0] = maxAmp2x2 ;   
318   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
319   triggerAmplitudes[2] = maxAmpnxn ;   
320   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
321   //printf("\n triggerAmplitudes ");
322   //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
323   //printf("\n");
324   //tr.Print("");
325   //
326   // Trigger jet staff
327   //
328   if(tr.GetNJetThreshold()>0) {
329     // Jet phi/eta
330     Int_t n0 = triggerPosition.GetSize();
331     const TH2F *hpatch = tr.GetJetMatrixE();
332     triggerPosition.Set(n0 + 2);
333     for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);   
334     // Add jet ampitudes
335     n0 = triggerAmplitudes.GetSize();
336     triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
337     Double_t *ampJet = tr.GetL1JetThresholds();
338     for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
339       triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
340     }
341   }
342   esd->AddEMCALTriggerPosition(triggerPosition);
343   esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
344   // Fill trigger hists
345   //  AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
346
347   //########################################
348   //##############Fill CaloCells###############
349   //########################################
350
351   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
352   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
353   if (!branchdig) { 
354     AliError("can't get the branch with the EMCAL digits !");
355     return;
356   }
357   branchdig->SetAddress(&digits);
358   digitsTree->GetEvent(0);
359   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
360   AliDebug(1,Form("%d digits",nDigits));
361
362   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
363   emcCells.CreateContainer(nDigits);
364   emcCells.SetType(AliESDCaloCells::kEMCALCell);
365   Float_t energy = 0;
366   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
367     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
368     if(dig->GetAmplitude() > 0 ){
369           energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
370           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
371                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
372                   idignew++;
373           }
374     }
375   }
376   emcCells.SetNumberOfCells(idignew);
377   emcCells.Sort();
378
379   //------------------------------------------------------------
380   //-----------------CLUSTERS-----------------------------
381   //------------------------------------------------------------
382   clustersTree->SetBranchStatus("*",0); //disable all branches
383   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
384
385   TObjArray *clusters = new TObjArray(100);
386   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
387   branch->SetAddress(&clusters);
388   branch->GetEntry(0);
389   //clustersTree->GetEvent(0);
390
391   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
392   AliDebug(1,Form("%d clusters",nClusters));
393
394   //######################################################
395   //#######################TRACK MATCHING###############
396   //######################################################
397   //Fill list of integers, each one is index of track to which the cluster belongs.
398
399   // step 1 - initialize array of matched track indexes
400   Int_t *matchedTrack = new Int_t[nClusters];
401   for (Int_t iclus = 0; iclus < nClusters; iclus++)
402     matchedTrack[iclus] = -1;  // neg. index --> no matched track
403   
404   // step 2, change the flag for all matched clusters found in tracks
405   Int_t iemcalMatch = -1;
406   Int_t endtpc = esd->GetNumberOfTracks();
407   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
408     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
409     iemcalMatch = track->GetEMCALcluster();
410     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
411   } 
412   
413   //########################################
414   //##############Fill CaloClusters#############
415   //########################################
416   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
417     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
418     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
419     if (Debug()) clust->Print();
420     // Get information from EMCAL reconstruction points
421     Float_t xyz[3];
422     TVector3 gpos;
423     clust->GetGlobalPosition(gpos);
424     for (Int_t ixyz=0; ixyz<3; ixyz++)
425       xyz[ixyz] = gpos[ixyz];
426     Float_t elipAxis[2];
427     clust->GetElipsAxis(elipAxis);
428        //Create digits lists
429     Int_t cellMult = clust->GetMultiplicity();
430     //TArrayS digiList(digitMult);
431     Float_t *amplFloat = clust->GetEnergiesList();
432     Int_t   *digitInts = clust->GetAbsId();
433     TArrayS absIdList(cellMult);
434     TArrayD fracList(cellMult);
435
436     Int_t newCellMult = 0;
437     for (Int_t iCell=0; iCell<cellMult; iCell++) {
438       if (amplFloat[iCell] > 0) {
439       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
440       //Uncomment when unfolding is done
441       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
442       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
443       //else
444       fracList[newCellMult] = 0; 
445       newCellMult++;
446       }
447     }
448
449     absIdList.Set(newCellMult);
450     fracList.Set(newCellMult);
451     
452     if(newCellMult > 0) { // accept cluster if it has some digit
453       nClustersNew++;
454       //Primaries
455       Int_t  parentMult  = 0;
456       Int_t *parentList =  clust->GetParents(parentMult);
457       // fills the ESDCaloCluster
458       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
459       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
460       ec->SetPosition(xyz);
461       ec->SetE(clust->GetEnergy());
462                 
463           //Distance to the nearest bad crystal
464           ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
465
466       ec->SetNCells(newCellMult);
467       //Change type of list from short to ushort
468       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
469       Double_t *newFracList  = new Double_t[newCellMult];
470       for(Int_t i = 0; i < newCellMult ; i++) {
471         newAbsIdList[i]=absIdList[i];
472         newFracList[i]=fracList[i];
473       }
474       ec->SetCellsAbsId(newAbsIdList);
475       ec->SetCellsAmplitudeFraction(newFracList);
476       ec->SetClusterDisp(clust->GetDispersion());
477       ec->SetClusterChi2(-1); //not yet implemented
478       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
479       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
480       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
481       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
482       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
483       arrayTrackMatched[0]= matchedTrack[iClust];
484       ec->AddTracksMatched(arrayTrackMatched);
485
486       TArrayI arrayParents(parentMult,parentList);
487       ec->AddLabels(arrayParents);
488
489       // add the cluster to the esd object
490       esd->AddCaloCluster(ec);
491       delete ec;
492       delete [] newAbsIdList ;
493       delete [] newFracList ;
494    }
495  } // cycle on clusters
496
497  delete [] matchedTrack;
498
499  //Fill ESDCaloCluster with PID weights
500  AliEMCALPID *pid = new AliEMCALPID;
501  //pid->SetPrintInfo(kTRUE);
502  pid->SetReconstructor(kTRUE);
503  pid->RunPID(esd);
504  delete pid;
505   
506  delete digits;
507  delete clusters;
508   
509  //Store EMCAL misalignment matrixes
510  FillMisalMatrixes(esd) ;
511
512 }
513
514 //==================================================================================
515 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
516         //Store EMCAL matrixes in ESD Header
517         
518         //Check, if matrixes was already stored
519         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
520                 if(esd->GetEMCALMatrix(sm)!=0)
521                         return ;
522         }
523         
524         //Create and store matrixes
525         if(!gGeoManager){
526                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
527                 return ;
528         }
529         //Note, that owner of copied marixes will be header
530         char path[255] ;
531         TGeoHMatrix * m = 0x0;
532         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
533                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
534                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
535                 
536                 if (gGeoManager->CheckPath(path)){
537                         gGeoManager->cd(path);
538                         m = gGeoManager->GetCurrentMatrix() ;
539 //                      printf("================================================= \n");
540 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
541 //                      m->Print("");
542                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
543 //                      printf("================================================= \n");
544                 }
545                 else{
546                         esd->SetEMCALMatrix(NULL,sm) ;
547                 }
548         }
549 }
550
551
552
553 //__________________________________________________________________________
554 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
555 {
556   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
557   if(fgDigitsArr) {
558     // Clear previous digits 
559     fgDigitsArr->Delete();
560     delete fgDigitsArr;
561   }
562   // Read the digits from the input tree
563   TBranch *branch = digitsTree->GetBranch("EMCAL");
564   if (!branch) { 
565     AliError("can't get the branch with the EMCAL digits !");
566     return;
567   }
568   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
569   branch->SetAddress(&fgDigitsArr);
570   branch->GetEntry(0);
571 }
572
573