]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Memory leak fixes:
[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, treeV0, 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         
241   fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
242
243   digitsTree->Fill();
244   digitsArr->Delete();
245   digitsTrg->Delete();
246   delete digitsArr;
247   delete digitsTrg;
248
249 }
250
251
252 //____________________________________________________________________________
253 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
254                                     AliESDEvent* esd) const
255 {
256   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
257   // and V0 
258   // Works on the current event
259   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
260   //return;
261
262   //######################################################
263   //#########Calculate trigger and set trigger info###########
264   //######################################################
265   // Obsolete, to be changed with new trigger emulator when consensus is achieved about what is stored in ESDs.
266   AliEMCALTrigger tr;
267   //   tr.SetPatchSize(1);  // create 4x4 patches
268   tr.SetSimulation(kFALSE); // Reconstruction mode
269   tr.SetDigitsList(fgDigitsArr);
270   // Get VZERO total multiplicity for jet trigger simulation 
271   // The simulation of jey trigger will be incorrect if no VZERO data 
272   // at ESD
273   AliESDVZERO* vZero = esd->GetVZEROData();
274   if(vZero) {
275     tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
276   }
277   //
278   tr.Trigger();
279
280   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
281   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
282   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
283   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
284
285   Int_t iSM2x2      = tr.Get2x2SuperModule();
286   Int_t iSMnxn      = tr.GetnxnSuperModule();
287   Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
288   Int_t iModulePhinxn = tr.GetnxnModulePhi();
289   Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
290   Int_t iModuleEtanxn = tr.GetnxnModuleEta();
291
292   AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d",  maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
293   AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d",  maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
294
295   TVector3    pos2x2(-1,-1,-1);
296   TVector3    posnxn(-1,-1,-1);
297
298   Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
299   Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
300   fGeom->GetGlobal(iAbsId2x2, pos2x2);
301   fGeom->GetGlobal(iAbsIdnxn, posnxn);
302   //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
303   
304   TArrayF triggerPosition(6);
305   triggerPosition[0] = pos2x2(0) ;   
306   triggerPosition[1] = pos2x2(1) ;   
307   triggerPosition[2] = pos2x2(2) ;  
308   triggerPosition[3] = posnxn(0) ;   
309   triggerPosition[4] = posnxn(1) ;   
310   triggerPosition[5] = posnxn(2) ;
311   //printf(" triggerPosition ");
312   //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
313
314   TArrayF triggerAmplitudes(4);
315   triggerAmplitudes[0] = maxAmp2x2 ;   
316   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
317   triggerAmplitudes[2] = maxAmpnxn ;   
318   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
319   //printf("\n triggerAmplitudes ");
320   //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
321   //printf("\n");
322   //tr.Print("");
323   //
324   // Trigger jet staff
325   //
326   if(tr.GetNJetThreshold()>0) {
327     // Jet phi/eta
328     Int_t n0 = triggerPosition.GetSize();
329     const TH2F *hpatch = tr.GetJetMatrixE();
330     triggerPosition.Set(n0 + 2);
331     for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);   
332     // Add jet ampitudes
333     n0 = triggerAmplitudes.GetSize();
334     triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
335     Double_t *ampJet = tr.GetL1JetThresholds();
336     for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
337       triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
338     }
339   }
340   esd->AddEMCALTriggerPosition(triggerPosition);
341   esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
342   // Fill trigger hists
343   //  AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
344
345   //########################################
346   //##############Fill CaloCells###############
347   //########################################
348
349   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
350   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
351   if (!branchdig) { 
352     AliError("can't get the branch with the EMCAL digits !");
353     return;
354   }
355   branchdig->SetAddress(&digits);
356   digitsTree->GetEvent(0);
357   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
358   AliDebug(1,Form("%d digits",nDigits));
359
360   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
361   emcCells.CreateContainer(nDigits);
362   emcCells.SetType(AliESDCaloCells::kEMCALCell);
363   Float_t energy = 0;
364   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
365     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
366     if(dig->GetAmplitude() > 0 ){
367           energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
368           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
369                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
370                   idignew++;
371           }
372     }
373   }
374   emcCells.SetNumberOfCells(idignew);
375   emcCells.Sort();
376
377   //------------------------------------------------------------
378   //-----------------CLUSTERS-----------------------------
379   //------------------------------------------------------------
380   clustersTree->SetBranchStatus("*",0); //disable all branches
381   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
382
383   TObjArray *clusters = new TObjArray(100);
384   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
385   branch->SetAddress(&clusters);
386   branch->GetEntry(0);
387   //clustersTree->GetEvent(0);
388
389   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
390   AliDebug(1,Form("%d clusters",nClusters));
391
392   //######################################################
393   //#######################TRACK MATCHING###############
394   //######################################################
395   //Fill list of integers, each one is index of track to which the cluster belongs.
396
397   // step 1 - initialize array of matched track indexes
398   Int_t *matchedTrack = new Int_t[nClusters];
399   for (Int_t iclus = 0; iclus < nClusters; iclus++)
400     matchedTrack[iclus] = -1;  // neg. index --> no matched track
401   
402   // step 2, change the flag for all matched clusters found in tracks
403   Int_t iemcalMatch = -1;
404   Int_t endtpc = esd->GetNumberOfTracks();
405   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
406     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
407     iemcalMatch = track->GetEMCALcluster();
408     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
409   } 
410   
411   //########################################
412   //##############Fill CaloClusters#############
413   //########################################
414   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
415     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
416     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
417     if (Debug()) clust->Print();
418     // Get information from EMCAL reconstruction points
419     Float_t xyz[3];
420     TVector3 gpos;
421     clust->GetGlobalPosition(gpos);
422     for (Int_t ixyz=0; ixyz<3; ixyz++)
423       xyz[ixyz] = gpos[ixyz];
424     Float_t elipAxis[2];
425     clust->GetElipsAxis(elipAxis);
426        //Create digits lists
427     Int_t cellMult = clust->GetMultiplicity();
428     //TArrayS digiList(digitMult);
429     Float_t *amplFloat = clust->GetEnergiesList();
430     Int_t   *digitInts = clust->GetAbsId();
431     TArrayS absIdList(cellMult);
432     TArrayD fracList(cellMult);
433
434     Int_t newCellMult = 0;
435     for (Int_t iCell=0; iCell<cellMult; iCell++) {
436       if (amplFloat[iCell] > 0) {
437       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
438       //Uncomment when unfolding is done
439       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
440       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
441       //else
442       fracList[newCellMult] = 0; 
443       newCellMult++;
444       }
445     }
446
447     absIdList.Set(newCellMult);
448     fracList.Set(newCellMult);
449     
450     if(newCellMult > 0) { // accept cluster if it has some digit
451       nClustersNew++;
452       //Primaries
453       Int_t  parentMult  = 0;
454       Int_t *parentList =  clust->GetParents(parentMult);
455       // fills the ESDCaloCluster
456       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
457       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
458       ec->SetPosition(xyz);
459       ec->SetE(clust->GetEnergy());
460                 
461           //Distance to the nearest bad crystal
462           ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
463
464       ec->SetNCells(newCellMult);
465       //Change type of list from short to ushort
466       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
467       Double_t *newFracList  = new Double_t[newCellMult];
468       for(Int_t i = 0; i < newCellMult ; i++) {
469         newAbsIdList[i]=absIdList[i];
470         newFracList[i]=fracList[i];
471       }
472       ec->SetCellsAbsId(newAbsIdList);
473       ec->SetCellsAmplitudeFraction(newFracList);
474       ec->SetClusterDisp(clust->GetDispersion());
475       ec->SetClusterChi2(-1); //not yet implemented
476       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
477       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
478       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
479       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
480       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
481       arrayTrackMatched[0]= matchedTrack[iClust];
482       ec->AddTracksMatched(arrayTrackMatched);
483
484       TArrayI arrayParents(parentMult,parentList);
485       ec->AddLabels(arrayParents);
486
487       // add the cluster to the esd object
488       esd->AddCaloCluster(ec);
489       delete ec;
490       delete [] newAbsIdList ;
491       delete [] newFracList ;
492    }
493  } // cycle on clusters
494
495  delete [] matchedTrack;
496
497  //Fill ESDCaloCluster with PID weights
498  AliEMCALPID *pid = new AliEMCALPID;
499  //pid->SetPrintInfo(kTRUE);
500  pid->SetReconstructor(kTRUE);
501  pid->RunPID(esd);
502  delete pid;
503   
504  delete digits;
505  delete clusters;
506   
507  //Store EMCAL misalignment matrixes
508  FillMisalMatrixes(esd) ;
509
510 }
511
512 //==================================================================================
513 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
514         //Store EMCAL matrixes in ESD Header
515         
516         //Check, if matrixes was already stored
517         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
518                 if(esd->GetEMCALMatrix(sm)!=0)
519                         return ;
520         }
521         
522         //Create and store matrixes
523         if(!gGeoManager){
524                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
525                 return ;
526         }
527         //Note, that owner of copied marixes will be header
528         char path[255] ;
529         TGeoHMatrix * m = 0x0;
530         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
531                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
532                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
533                 
534                 if (gGeoManager->CheckPath(path)){
535                         gGeoManager->cd(path);
536                         m = gGeoManager->GetCurrentMatrix() ;
537 //                      printf("================================================= \n");
538 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
539 //                      m->Print("");
540                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
541 //                      printf("================================================= \n");
542                 }
543                 else{
544                         esd->SetEMCALMatrix(NULL,sm) ;
545                 }
546         }
547 }
548
549
550
551 //__________________________________________________________________________
552 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
553 {
554   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
555   if(fgDigitsArr) {
556     // Clear previous digits 
557     fgDigitsArr->Delete();
558     delete fgDigitsArr;
559   }
560   // Read the digits from the input tree
561   TBranch *branch = digitsTree->GetBranch("EMCAL");
562   if (!branch) { 
563     AliError("can't get the branch with the EMCAL digits !");
564     return;
565   }
566   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
567   branch->SetAddress(&fgDigitsArr);
568   branch->GetEntry(0);
569 }
570
571