]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
bug in loop, wrong varialble in increment
[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 "AliCaloCalibPedestal.h"
41 #include "AliEMCALCalibData.h"
42 #include "AliESDEvent.h"
43 #include "AliESDCaloCluster.h"
44 #include "AliESDCaloCells.h"
45 #include "AliESDtrack.h"
46 #include "AliEMCALLoader.h"
47 #include "AliEMCALRawUtils.h"
48 #include "AliEMCALDigit.h"
49 #include "AliEMCALClusterizerv1.h"
50 #include "AliEMCALClusterizerNxN.h"
51 #include "AliEMCALRecPoint.h"
52 #include "AliEMCALPID.h"
53 #include "AliEMCALTrigger.h"
54 #include "AliRawReader.h"
55 #include "AliCDBEntry.h"
56 #include "AliCDBManager.h"
57 #include "AliEMCALGeometry.h"
58 #include "AliEMCAL.h"
59 #include "AliESDVZERO.h"
60 #include "AliCDBManager.h"
61 #include "AliRunLoader.h"
62 #include "AliRun.h"
63 #include "AliEMCALTriggerData.h"
64 #include "AliEMCALTriggerElectronics.h"
65 #include "AliEMCALTriggerDCSConfigDB.h"
66 #include "AliEMCALTriggerDCSConfig.h"
67 #include "AliEMCALTriggerData.h"
68 #include "AliEMCALTriggerRawDigit.h"
69 #include "AliEMCALTriggerPatch.h"
70 #include "AliEMCALTriggerTypes.h"
71
72 ClassImp(AliEMCALReconstructor) 
73
74 const AliEMCALRecParam*     AliEMCALReconstructor::fgkRecParam        = 0;   // EMCAL rec. parameters
75 AliEMCALRawUtils*           AliEMCALReconstructor::fgRawUtils         = 0;   // EMCAL raw utilities class
76 AliEMCALClusterizer*        AliEMCALReconstructor::fgClusterizer      = 0;   // EMCAL clusterizer class
77 TClonesArray*               AliEMCALReconstructor::fgDigitsArr        = 0;   // list of digits, to be used multiple times
78 TObjArray*                  AliEMCALReconstructor::fgClustersArr      = 0;   // list of clusters, to be used multiple times
79 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
80 //____________________________________________________________________________
81 AliEMCALReconstructor::AliEMCALReconstructor() 
82   : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0),fTriggerData(0x0) 
83 {
84   // ctor
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   InitClusterizer();
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         
140
141
142 //____________________________________________________________________________
143 AliEMCALReconstructor::~AliEMCALReconstructor()
144 {
145   // dtor
146
147   if(fGeom)              delete fGeom;
148   if(fCalibData)         delete fCalibData;
149   if(fPedestalData)      delete fPedestalData;
150   
151   if(fgDigitsArr){
152     fgDigitsArr->Clear("C");
153     delete fgDigitsArr; 
154   }
155   
156   if(fgClustersArr){
157     fgClustersArr->Clear();
158     delete fgClustersArr; 
159   }
160   
161   if(fgRawUtils)         delete fgRawUtils;
162   if(fgClusterizer)      delete fgClusterizer;
163   if(fgTriggerProcessor) delete fgTriggerProcessor;
164   
165   AliCodeTimer::Instance()->Print();
166
167
168 // //____________________________________________________________________________
169 // void AliEMCALReconstructor::Init()
170 // {
171 //   // Trigger hists - Oct 24, 2007
172 //    fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
173 // }
174
175 //____________________________________________________________________________
176 void AliEMCALReconstructor::InitClusterizer() 
177 {
178   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
179   
180   AliEMCALRecParam *recParam = NULL;
181   AliCDBEntry *entry = (AliCDBEntry*) 
182   AliCDBManager::Instance()->Get("EMCAL/Calib/RecoParam");
183   //Get The reco param for the default event specie
184   if (entry) 
185     recParam = (AliEMCALRecParam*)((TObjArray *) entry->GetObject())->At(0);
186   
187   if(!recParam){  
188     AliFatal("RecoParam not found in CDB!");
189   }
190   else{
191     if (recParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
192     {
193       fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData); 
194     }
195     else
196     {
197       fgClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData,fPedestalData); 
198     }
199   }
200   
201 }
202
203 //____________________________________________________________________________
204 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
205 {
206   // method called by AliReconstruction; 
207   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
208   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
209   // the global tracking.
210   // Works on the current event.
211
212   AliCodeTimerAuto("",0)
213
214   ReadDigitsArrayFromTree(digitsTree);
215
216   fgClusterizer->InitParameters();
217   fgClusterizer->SetOutput(clustersTree);
218         
219   //Skip clusterization of LED events
220   if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
221
222                 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
223
224                   fgClusterizer->SetInput(digitsTree);
225     
226                   if(Debug())
227                           fgClusterizer->Digits2Clusters("deb all") ;
228                   else
229                           fgClusterizer->Digits2Clusters("");
230     
231                   fgClusterizer->Clear();
232
233           }//digits array exists and has somethind
234   }//not a LED event
235         
236   clustersTree->Fill(); 
237 }
238
239 //____________________________________________________________________________
240 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
241
242 {
243   // Conversion from raw data to
244   // EMCAL digits.
245   // Works on a single-event basis
246
247   rawReader->Reset() ; 
248
249   fTriggerData->SetMode(1);     
250
251   if(fgDigitsArr) fgDigitsArr->Clear("C");
252
253   TClonesArray *digitsTrg = new TClonesArray("AliEMCALTriggerRawDigit", 32 * 96);
254
255   Int_t bufsize = 32000;
256   digitsTree->Branch("EMCAL", &fgDigitsArr, bufsize);
257   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
258         
259   //Skip calibration events do the rest
260   Bool_t doFit = kTRUE;
261   if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
262   if (doFit){
263           //must be done here because, in constructor, option is not yet known
264           fgRawUtils->SetOption(GetOption());
265
266           fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
267           fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
268           fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
269           fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
270           fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
271           fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
272           fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
273           fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
274           fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
275           fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
276         
277           fgRawUtils->Raw2Digits(rawReader,fgDigitsArr,fPedestalData,digitsTrg,fTriggerData);
278   }//skip calibration event
279   else{
280         AliDebug(1," Calibration Event, skip!");
281   }
282         
283   digitsTree->Fill();
284   digitsTrg->Delete();
285   delete digitsTrg;
286
287 }
288
289
290 //____________________________________________________________________________
291 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
292                                     AliESDEvent* esd) const
293 {
294   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
295   // and V0 
296   // Works on the current event
297   // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
298   //return;
299
300   //FIXME UNCOMMENT WHEN ESDTRIGGER AVAILABLE 
301 //   // Trigger
302 //   Int_t v0M[2] = {0, 0};
303         
304 //   AliESDVZERO* esdV0 = esd->GetVZEROData();
305
306 //   if (esdV0) 
307 //   {
308 //        for (Int_t i = 0; i < 32; i++)
309 //        {
310 //                v0M[0] += esdV0->GetAdcV0C(i);
311 //                v0M[1] += esdV0->GetAdcV0A(i);
312 //        }
313 //   }
314 //   else
315 //   {
316 //        AliWarning("Cannot retrieve V0 ESD! Run w/ null V0 charges");
317 //   }
318
319 //   TClonesArray *trgDigits = new TClonesArray("AliEMCALTriggerRawDigit",1000);
320         
321 //   TBranch *branchtrg = digitsTree->GetBranch("EMTRG");
322         
323 //   if (!branchtrg) 
324 //   { 
325 //        AliError("Can't get the branch with the EMCAL trigger digits!");
326 //        return;
327 //   }
328         
329 //   branchtrg->SetAddress(&trgDigits);
330 //   branchtrg->GetEntry(0);
331   
332 //   // Note: fgTriggerProcessor reset done at the end of this method
333 //   fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, fTriggerData);
334
335 //   // Fill ESD
336 //   AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("EMCAL");
337   
338 //   if (trgESD)
339 //   {
340 //        trgESD->Allocate(trgDigits->GetEntriesFast());
341           
342 //        for (Int_t i = 0; i < trgDigits->GetEntriesFast(); i++)
343 //        {       
344 //                AliEMCALTriggerRawDigit* rdig = (AliEMCALTriggerRawDigit*)trgDigits->At(i);
345                   
346 //                Int_t px, py;
347 //                if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
348 //                {
349 //                        Int_t a = -1, t = -1, times[10]; 
350                           
351 //                        rdig->GetMaximum(a, t);
352 //                        rdig->GetL0Times(times);
353                           
354 // //                     rdig->Print("");
355                           
356 //                        trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum());
357 //                }
358 //        }
359           
360 // //     cout << "End of Adding................." << endl;
361
362 //        trgESD->SetL1Threshold(0, fTriggerData->GetL1GammaThreshold());
363           
364 //        trgESD->SetL1Threshold(1, fTriggerData->GetL1JetThreshold()  );
365
366 //        for (Int_t i = 0; i < kTriggerTypeEnd; i++)
367 //        {       
368 //                for (Int_t j = 0; j < 2; j++)
369 //                {
370 //                        TClonesArray* patches = fTriggerData->GetPatches((TriggerType_t)i, j);
371                           
372 //                        TIter NextPatch(patches);
373 //                        while (AliEMCALTriggerPatch* p = (AliEMCALTriggerPatch*)NextPatch())
374 //                        {
375 //                                TVector2 pos; p->Position(pos);
376 //                                trgESD->SetTriggerBits(pos.X(), pos.Y(), i, j);
377 //                        }
378 //                }
379 //        }
380 //   }
381
382 //   // Resetting
383 //   fTriggerData->Reset();
384 // //  cout << "Reset trg data" << endl;
385   //FIXME UNCOMMENT WHEN ESDTRIGGER AVAILABLE 
386   
387   //########################################
388   //##############Fill CaloCells###############
389   //########################################
390   ReadDigitsArrayFromTree(digitsTree);
391
392 //   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
393 //   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
394 //   if (!branchdig) { 
395 //     AliError("can't get the branch with the EMCAL digits !");
396 //     return;
397 //   }
398 //   branchdig->SetAddress(&digits);
399 //   digitsTree->GetEvent(0);
400   Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
401   AliDebug(1,Form("%d digits",nDigits));
402
403   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
404   emcCells.CreateContainer(nDigits);
405   emcCells.SetType(AliVCaloCells::kEMCALCell);
406   Float_t energy = 0;
407   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
408     const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
409     if(dig->GetAmplitude() > 0 ){
410           energy = fgClusterizer->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
411           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
412                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
413                   idignew++;
414           }
415     }
416   }
417   emcCells.SetNumberOfCells(idignew);
418   emcCells.Sort();
419
420   //------------------------------------------------------------
421   //-----------------CLUSTERS-----------------------------
422   //------------------------------------------------------------
423   clustersTree->SetBranchStatus("*",0); //disable all branches
424   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
425   if(fgClustersArr) fgClustersArr->Clear();
426   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
427   branch->SetAddress(&fgClustersArr);
428   branch->GetEntry(0);
429   //clustersTree->GetEvent(0);
430
431   Int_t nClusters = fgClustersArr->GetEntries(),  nClustersNew=0;
432   AliDebug(1,Form("%d clusters",nClusters));
433
434   //######################################################
435   //#######################TRACK MATCHING###############
436   //######################################################
437   //Fill list of integers, each one is index of track to which the cluster belongs.
438
439   // step 1 - initialize array of matched track indexes
440   Int_t *matchedTrack = new Int_t[nClusters];
441   for (Int_t iclus = 0; iclus < nClusters; iclus++)
442     matchedTrack[iclus] = -1;  // neg. index --> no matched track
443   
444   // step 2, change the flag for all matched clusters found in tracks
445   Int_t iemcalMatch = -1;
446   Int_t endtpc = esd->GetNumberOfTracks();
447   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
448     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
449     iemcalMatch = track->GetEMCALcluster();
450     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
451   } 
452   
453   //########################################
454   //##############Fill CaloClusters#############
455   //########################################
456   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
457     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
458     //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
459     if (Debug()) clust->Print();
460     // Get information from EMCAL reconstruction points
461     Float_t xyz[3];
462     TVector3 gpos;
463     clust->GetGlobalPosition(gpos);
464     for (Int_t ixyz=0; ixyz<3; ixyz++)
465       xyz[ixyz] = gpos[ixyz];
466     Float_t elipAxis[2];
467     clust->GetElipsAxis(elipAxis);
468        //Create digits lists
469     Int_t cellMult = clust->GetMultiplicity();
470     //TArrayS digiList(digitMult);
471     Float_t *amplFloat = clust->GetEnergiesList();
472     Int_t   *digitInts = clust->GetAbsId();
473     TArrayS absIdList(cellMult);
474     TArrayD fracList(cellMult);
475
476     Int_t newCellMult = 0;
477     for (Int_t iCell=0; iCell<cellMult; iCell++) {
478       if (amplFloat[iCell] > 0) {
479       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
480       //Uncomment when unfolding is done
481       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
482       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
483       //else
484       fracList[newCellMult] = 0; 
485       newCellMult++;
486       }
487     }
488
489     absIdList.Set(newCellMult);
490     fracList.Set(newCellMult);
491     
492     if(newCellMult > 0) { // accept cluster if it has some digit
493       nClustersNew++;
494       //Primaries
495       Int_t  parentMult  = 0;
496       Int_t *parentList =  clust->GetParents(parentMult);
497       // fills the ESDCaloCluster
498       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
499       ec->SetType(AliVCluster::kEMCALClusterv1);
500       ec->SetPosition(xyz);
501       ec->SetE(clust->GetEnergy());
502                 
503       //Distance to the nearest bad crystal
504       ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
505
506       ec->SetNCells(newCellMult);
507       //Change type of list from short to ushort
508       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
509       Double_t *newFracList   = new Double_t[newCellMult];
510       for(Int_t i = 0; i < newCellMult ; i++) {
511         newAbsIdList[i]=absIdList[i];
512         newFracList[i] =fracList[i];
513       }
514       ec->SetCellsAbsId(newAbsIdList);
515       ec->SetCellsAmplitudeFraction(newFracList);
516       ec->SetDispersion(clust->GetDispersion());
517       ec->SetChi2(-1); //not yet implemented
518       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
519       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
520       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
521       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
522       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
523       arrayTrackMatched[0]= matchedTrack[iClust];
524       ec->AddTracksMatched(arrayTrackMatched);
525
526       TArrayI arrayParents(parentMult,parentList);
527       ec->AddLabels(arrayParents);
528
529       // add the cluster to the esd object
530       esd->AddCaloCluster(ec);
531       delete ec;
532       delete [] newAbsIdList ;
533       delete [] newFracList ;
534    }
535  } // cycle on clusters
536
537  delete [] matchedTrack;
538
539  //Fill ESDCaloCluster with PID weights
540  AliEMCALPID *pid = new AliEMCALPID;
541  //pid->SetPrintInfo(kTRUE);
542  pid->SetReconstructor(kTRUE);
543  pid->RunPID(esd);
544  delete pid;
545     
546  //Store EMCAL misalignment matrixes
547  FillMisalMatrixes(esd) ;
548
549 }
550
551 //==================================================================================
552 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
553         //Store EMCAL matrixes in ESD Header
554         
555         //Check, if matrixes was already stored
556         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
557                 if(esd->GetEMCALMatrix(sm)!=0)
558                         return ;
559         }
560         
561         //Create and store matrixes
562         if(!gGeoManager){
563                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
564                 return ;
565         }
566         //Note, that owner of copied marixes will be header
567   const Int_t bufsize = 255;
568         char path[bufsize] ;
569         TGeoHMatrix * m = 0x0;
570         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
571                 snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
572                 if(sm >= 10) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
573                 
574                 if (gGeoManager->CheckPath(path)){
575                         gGeoManager->cd(path);
576                         m = gGeoManager->GetCurrentMatrix() ;
577 //                      printf("================================================= \n");
578 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
579 //                      m->Print("");
580                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
581 //                      printf("================================================= \n");
582                 }
583                 else{
584                         esd->SetEMCALMatrix(NULL,sm) ;
585                 }
586         }
587 }
588
589
590
591 //__________________________________________________________________________
592 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
593 {
594   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
595   if(fgDigitsArr) {
596     // Clear previous digits 
597     fgDigitsArr->Clear("C");
598     //delete fgDigitsArr;
599   }
600   // Read the digits from the input tree
601   TBranch *branch = digitsTree->GetBranch("EMCAL");
602   if (!branch) { 
603     AliError("can't get the branch with the EMCAL digits !");
604     return;
605   }
606   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
607   branch->SetAddress(&fgDigitsArr);
608   branch->GetEntry(0);
609 }
610
611