]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Minor 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 "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   //########################################
303   // Trigger
304   //########################################
305
306    Int_t v0M[2] = {0, 0};
307         
308    AliESDVZERO* esdV0 = esd->GetVZEROData();
309
310    if (esdV0) 
311    {
312           for (Int_t i = 0; i < 32; i++)
313           {
314                   v0M[0] += (Int_t)esdV0->GetAdcV0C(i);
315                   v0M[1] += (Int_t)esdV0->GetAdcV0A(i);
316           }
317    }
318    else
319    {
320           AliWarning("Cannot retrieve V0 ESD! Run w/ null V0 charges");
321    }
322
323    TClonesArray *trgDigits = new TClonesArray("AliEMCALTriggerRawDigit",1000);
324         
325    TBranch *branchtrg = digitsTree->GetBranch("EMTRG");
326         
327    if (!branchtrg) 
328    { 
329           AliError("Can't get the branch with the EMCAL trigger digits!");
330           return;
331    }
332         
333    branchtrg->SetAddress(&trgDigits);
334    branchtrg->GetEntry(0);
335   
336    // Note: fgTriggerProcessor reset done at the end of this method
337    fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, fTriggerData);
338
339    // Fill ESD
340    AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("EMCAL");
341   
342    if (trgESD)
343    {
344           trgESD->Allocate(trgDigits->GetEntriesFast());
345           
346           for (Int_t i = 0; i < trgDigits->GetEntriesFast(); i++)
347           {       
348                   AliEMCALTriggerRawDigit* rdig = (AliEMCALTriggerRawDigit*)trgDigits->At(i);
349                   
350                   Int_t px, py;
351                   if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
352                   {
353                           Int_t a = -1, t = -1, times[10]; 
354                           
355                           rdig->GetMaximum(a, t);
356                           rdig->GetL0Times(times);
357                                                   
358                           trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum());
359                   }
360           }
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((Int_t) pos.X(), (Int_t) pos.Y(), i, j);
377                           }
378                   }
379           }
380    }
381
382    // Resetting
383    fTriggerData->Reset();
384   
385   //########################################
386   //##############Fill CaloCells###############
387   //########################################
388   ReadDigitsArrayFromTree(digitsTree);
389
390 //   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
391 //   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
392 //   if (!branchdig) { 
393 //     AliError("can't get the branch with the EMCAL digits !");
394 //     return;
395 //   }
396 //   branchdig->SetAddress(&digits);
397 //   digitsTree->GetEvent(0);
398   Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
399   AliDebug(1,Form("%d digits",nDigits));
400
401   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
402   emcCells.CreateContainer(nDigits);
403   emcCells.SetType(AliVCaloCells::kEMCALCell);
404   Float_t energy = 0;
405   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
406     const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
407     if(dig->GetAmplitude() > 0 ){
408           energy = fgClusterizer->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
409           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
410                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
411                   idignew++;
412           }
413     }
414   }
415   emcCells.SetNumberOfCells(idignew);
416   emcCells.Sort();
417
418   //------------------------------------------------------------
419   //-----------------CLUSTERS-----------------------------
420   //------------------------------------------------------------
421   clustersTree->SetBranchStatus("*",0); //disable all branches
422   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
423   if(fgClustersArr) fgClustersArr->Clear();
424   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
425   branch->SetAddress(&fgClustersArr);
426   branch->GetEntry(0);
427   //clustersTree->GetEvent(0);
428
429   Int_t nClusters = fgClustersArr->GetEntries(),  nClustersNew=0;
430   AliDebug(1,Form("%d clusters",nClusters));
431
432   //######################################################
433   //#######################TRACK MATCHING###############
434   //######################################################
435   //Fill list of integers, each one is index of track to which the cluster belongs.
436
437   // step 1 - initialize array of matched track indexes
438   Int_t *matchedTrack = new Int_t[nClusters];
439   for (Int_t iclus = 0; iclus < nClusters; iclus++)
440     matchedTrack[iclus] = -1;  // neg. index --> no matched track
441   
442   // step 2, change the flag for all matched clusters found in tracks
443   Int_t iemcalMatch = -1;
444   Int_t endtpc = esd->GetNumberOfTracks();
445   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
446     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
447     iemcalMatch = track->GetEMCALcluster();
448     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
449   } 
450   
451   //########################################
452   //##############Fill CaloClusters#############
453   //########################################
454   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
455     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
456     //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
457     if (Debug()) clust->Print();
458     // Get information from EMCAL reconstruction points
459     Float_t xyz[3];
460     TVector3 gpos;
461     clust->GetGlobalPosition(gpos);
462     for (Int_t ixyz=0; ixyz<3; ixyz++)
463       xyz[ixyz] = gpos[ixyz];
464     Float_t elipAxis[2];
465     clust->GetElipsAxis(elipAxis);
466        //Create digits lists
467     Int_t cellMult = clust->GetMultiplicity();
468     //TArrayS digiList(digitMult);
469     Float_t *amplFloat = clust->GetEnergiesList();
470     Int_t   *digitInts = clust->GetAbsId();
471     TArrayS absIdList(cellMult);
472     TArrayD fracList(cellMult);
473
474     Int_t newCellMult = 0;
475     for (Int_t iCell=0; iCell<cellMult; iCell++) {
476       if (amplFloat[iCell] > 0) {
477       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
478       //Calculate Fraction
479       if(emcCells.GetCellAmplitude(digitInts[iCell])>0 && GetRecParam()->GetUnfold())
480         fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell]));//get cell calibration value 
481       else
482         fracList[newCellMult] = 0; 
483       newCellMult++;
484       }
485     }
486
487     absIdList.Set(newCellMult);
488     fracList.Set(newCellMult);
489     
490     if(newCellMult > 0) { // accept cluster if it has some digit
491       nClustersNew++;
492       //Primaries
493       Int_t  parentMult  = 0;
494       Int_t *parentList =  clust->GetParents(parentMult);
495       // fills the ESDCaloCluster
496       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
497       ec->SetType(AliVCluster::kEMCALClusterv1);
498       ec->SetPosition(xyz);
499       ec->SetE(clust->GetEnergy());
500                 
501       //Distance to the nearest bad crystal
502       ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
503
504       ec->SetNCells(newCellMult);
505       //Change type of list from short to ushort
506       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
507       Double_t *newFracList   = new Double_t[newCellMult];
508       for(Int_t i = 0; i < newCellMult ; i++) {
509         newAbsIdList[i]=absIdList[i];
510         newFracList[i] =fracList[i];
511       }
512       ec->SetCellsAbsId(newAbsIdList);
513       ec->SetCellsAmplitudeFraction(newFracList);
514       ec->SetDispersion(clust->GetDispersion());
515       ec->SetChi2(-1); //not yet implemented
516       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
517       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
518       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
519       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
520       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
521       arrayTrackMatched[0]= matchedTrack[iClust];
522       ec->AddTracksMatched(arrayTrackMatched);
523
524       TArrayI arrayParents(parentMult,parentList);
525       ec->AddLabels(arrayParents);
526
527       // add the cluster to the esd object
528       esd->AddCaloCluster(ec);
529       delete ec;
530       delete [] newAbsIdList ;
531       delete [] newFracList ;
532    }
533  } // cycle on clusters
534
535  delete [] matchedTrack;
536
537  //Fill ESDCaloCluster with PID weights
538  AliEMCALPID *pid = new AliEMCALPID;
539  //pid->SetPrintInfo(kTRUE);
540  pid->SetReconstructor(kTRUE);
541  pid->RunPID(esd);
542  delete pid;
543     
544  //Store EMCAL misalignment matrixes
545  FillMisalMatrixes(esd) ;
546
547 }
548
549 //==================================================================================
550 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
551         //Store EMCAL matrixes in ESD Header
552         
553         //Check, if matrixes was already stored
554         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
555                 if(esd->GetEMCALMatrix(sm)!=0)
556                         return ;
557         }
558         
559         //Create and store matrixes
560         if(!gGeoManager){
561                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
562                 return ;
563         }
564         //Note, that owner of copied marixes will be header
565   const Int_t bufsize = 255;
566         char path[bufsize] ;
567         TGeoHMatrix * m = 0x0;
568         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
569                 snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
570                 if(sm >= 10) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
571                 
572                 if (gGeoManager->CheckPath(path)){
573                         gGeoManager->cd(path);
574                         m = gGeoManager->GetCurrentMatrix() ;
575 //                      printf("================================================= \n");
576 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
577 //                      m->Print("");
578                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
579 //                      printf("================================================= \n");
580                 }
581                 else{
582                         esd->SetEMCALMatrix(NULL,sm) ;
583                 }
584         }
585 }
586
587
588
589 //__________________________________________________________________________
590 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
591 {
592   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
593   if(fgDigitsArr) {
594     // Clear previous digits 
595     fgDigitsArr->Clear("C");
596     //delete fgDigitsArr;
597   }
598   // Read the digits from the input tree
599   TBranch *branch = digitsTree->GetBranch("EMCAL");
600   if (!branch) { 
601     AliError("can't get the branch with the EMCAL digits !");
602     return;
603   }
604   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
605   branch->SetAddress(&fgDigitsArr);
606   branch->GetEntry(0);
607 }
608
609