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