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