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