added separate histogram entry for V0A and C BG flag
[u/mrichter/AliRoot.git] / EMCAL / AliEMCALReconstructor.cxx
CommitLineData
f6019cda 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//_________________________________________________________________________
fa42b1f3 19//--
20//-- Yves Schutz (SUBATECH)
f6019cda 21// Reconstruction class. Redesigned from the old AliReconstructionner class and
22// derived from STEER/AliReconstructor.
23//
85c25c2e 24//-- Aleksei Pavlinov : added staf for EMCAL jet trigger 9Apr 25, 2008)
25// : fgDigitsArr should read just once at event
26
f6019cda 27// --- ROOT system ---
85c25c2e 28#include <TList.h>
29#include <TClonesArray.h>
60d8ad94 30#include <TH2.h>
0c5b726e 31#include "TGeoManager.h"
32#include "TGeoMatrix.h"
f6019cda 33
34// --- Standard library ---
35
36// --- AliRoot header files ---
f6019cda 37#include "AliEMCALReconstructor.h"
5dee926e 38
aaa3cb7c 39#include "AliCodeTimer.h"
af885e0f 40#include "AliESDEvent.h"
89ffc0b0 41#include "AliESDCaloCluster.h"
0e7c6655 42#include "AliESDCaloCells.h"
6a0cf740 43#include "AliESDtrack.h"
5dee926e 44#include "AliEMCALLoader.h"
98e9578e 45#include "AliEMCALRawUtils.h"
0e7c6655 46#include "AliEMCALDigit.h"
f6019cda 47#include "AliEMCALClusterizerv1.h"
5dee926e 48#include "AliEMCALRecPoint.h"
dc293ae9 49#include "AliEMCALPID.h"
0964c2e9 50#include "AliEMCALTrigger.h"
1d59832c 51#include "AliRawReader.h"
fa42b1f3 52#include "AliCDBEntry.h"
53#include "AliCDBManager.h"
65bdc82f 54#include "AliEMCALGeometry.h"
ac8ae9fe 55#include "AliEMCAL.h"
85c25c2e 56#include "AliEMCALHistoUtilities.h"
57#include "AliESDVZERO.h"
0c5b726e 58#include "AliCDBManager.h"
72c58de0 59#include "AliRunLoader.h"
60#include "AliRun.h"
1d59832c 61
85c25c2e 62ClassImp(AliEMCALReconstructor)
f6019cda 63
ba6de5ea 64const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
65bdc82f 65AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
9517d886 66AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0; // EMCAL clusterizer class
85c25c2e 67TClonesArray* AliEMCALReconstructor::fgDigitsArr = 0; // shoud read just once at event
f6019cda 68//____________________________________________________________________________
18a21c7c 69AliEMCALReconstructor::AliEMCALReconstructor()
40164976 70 : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0)
f6019cda 71{
72 // ctor
65bdc82f 73
74 fgRawUtils = new AliEMCALRawUtils;
72c58de0 75
76 //To make sure we match with the geometry in a simulation file,
77 //let's try to get it first. If not, take the default geometry
33c3c91a 78 AliRunLoader *rl = AliRunLoader::Instance();
72c58de0 79 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
80 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
81 } else {
82 AliInfo(Form("Using default geometry in reconstruction"));
937d0661 83 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
65bdc82f 84 }
0e7c6655 85
0c5b726e 86 //Get calibration parameters
87 if(!fCalibData)
88 {
89 AliCDBEntry *entry = (AliCDBEntry*)
90 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
91 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
92 }
93
94 if(!fCalibData)
95 AliFatal("Calibration parameters not found in CDB!");
96
40164976 97 //Get calibration parameters
98 if(!fPedestalData)
99 {
100 AliCDBEntry *entry = (AliCDBEntry*)
101 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
102 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
103 }
104
105 if(!fPedestalData)
106 AliFatal("Dead map not found in CDB!");
107
108
0c5b726e 109 //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
40164976 110 fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData);
0c5b726e 111
72c58de0 112 if(!fGeom) AliFatal(Form("Could not get geometry!"));
113
f6019cda 114}
115
0a4cb131 116//____________________________________________________________________________
f6019cda 117AliEMCALReconstructor::~AliEMCALReconstructor()
118{
119 // dtor
65bdc82f 120 delete fGeom;
aaa3cb7c 121 AliCodeTimer::Instance()->Print();
f6019cda 122}
123
124//____________________________________________________________________________
85c25c2e 125void AliEMCALReconstructor::Init()
126{
127 // Trigger hists - Oct 24, 2007
128 fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
129}
130
131//____________________________________________________________________________
c47157cd 132void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
f6019cda 133{
134 // method called by AliReconstruction;
135 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
136 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
137 // the global tracking.
c47157cd 138 // Works on the current event.
fa42b1f3 139
1a429a6b 140 AliCodeTimerAuto("",0)
aaa3cb7c 141
85c25c2e 142 ReadDigitsArrayFromTree(digitsTree);
0832a2bf 143 fgClusterizer->InitParameters();
4601e3a7 144 fgClusterizer->SetOutput(clustersTree);
0c5b726e 145
85c25c2e 146 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
147
9517d886 148 fgClusterizer->SetInput(digitsTree);
9517d886 149
85c25c2e 150 if(Debug())
9517d886 151 fgClusterizer->Digits2Clusters("deb all") ;
85c25c2e 152 else
9517d886 153 fgClusterizer->Digits2Clusters("");
154
155 fgClusterizer->Clear();
85c25c2e 156
157 }
9517d886 158
f6019cda 159}
160
161//____________________________________________________________________________
c47157cd 162void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
98e9578e 163
a68156e6 164{
c47157cd 165 // Conversion from raw data to
166 // EMCAL digits.
167 // Works on a single-event basis
85c60a8e 168
98e9578e 169 rawReader->Reset() ;
98e9578e 170
3e5fa09c 171 TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
c47157cd 172 Int_t bufsize = 32000;
173 digitsTree->Branch("EMCAL", &digitsArr, bufsize);
98e9578e 174
65bdc82f 175 //must be done here because, in constructor, option is not yet known
176 fgRawUtils->SetOption(GetOption());
b4133f05 177
0832a2bf 178 fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
179 fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
180 fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
181 fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
182 fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
b4133f05 183
65bdc82f 184 fgRawUtils->Raw2Digits(rawReader,digitsArr);
c615db53 185
186 digitsTree->Fill();
187 digitsArr->Delete();
188 delete digitsArr;
189
a68156e6 190}
191
85c25c2e 192
a68156e6 193//____________________________________________________________________________
0e7c6655 194void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
c47157cd 195 AliESDEvent* esd) const
f6019cda 196{
98e9578e 197 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
85c25c2e 198 // and V0
c47157cd 199 // Works on the current event
85c25c2e 200 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
201 //return;
92da3372 202
6a0cf740 203 //######################################################
0964c2e9 204 //#########Calculate trigger and set trigger info###########
6a0cf740 205 //######################################################
85c25c2e 206
207 AliEMCALTrigger tr;
208 // tr.SetPatchSize(1); // create 4x4 patches
209 tr.SetSimulation(kFALSE); // Reconstruction mode
210 tr.SetDigitsList(fgDigitsArr);
211 // Get VZERO total multiplicity for jet trigger simulation
212 // The simulation of jey trigger will be incorrect if no VZERO data
213 // at ESD
214 AliESDVZERO* vZero = esd->GetVZEROData();
215 if(vZero) {
216 tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
217 }
218 //
0964c2e9 219 tr.Trigger();
85c25c2e 220
0964c2e9 221 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
222 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
223 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
224 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
225
0964c2e9 226 Int_t iSM2x2 = tr.Get2x2SuperModule();
227 Int_t iSMnxn = tr.GetnxnSuperModule();
85c25c2e 228 Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
229 Int_t iModulePhinxn = tr.GetnxnModulePhi();
230 Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
231 Int_t iModuleEtanxn = tr.GetnxnModuleEta();
0964c2e9 232
85c25c2e 233 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
234 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
0964c2e9 235
236 TVector3 pos2x2(-1,-1,-1);
237 TVector3 posnxn(-1,-1,-1);
238
85c25c2e 239 Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
240 Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
65bdc82f 241 fGeom->GetGlobal(iAbsId2x2, pos2x2);
242 fGeom->GetGlobal(iAbsIdnxn, posnxn);
85c25c2e 243 //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
0964c2e9 244
245 TArrayF triggerPosition(6);
246 triggerPosition[0] = pos2x2(0) ;
247 triggerPosition[1] = pos2x2(1) ;
248 triggerPosition[2] = pos2x2(2) ;
249 triggerPosition[3] = posnxn(0) ;
250 triggerPosition[4] = posnxn(1) ;
85c25c2e 251 triggerPosition[5] = posnxn(2) ;
252 //printf(" triggerPosition ");
253 //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
0964c2e9 254
255 TArrayF triggerAmplitudes(4);
256 triggerAmplitudes[0] = maxAmp2x2 ;
257 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
258 triggerAmplitudes[2] = maxAmpnxn ;
259 triggerAmplitudes[3] = ampOutOfPatchnxn ;
85c25c2e 260 //printf("\n triggerAmplitudes ");
261 //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
262 //printf("\n");
85d4cbde 263 //tr.Print("");
60d8ad94 264 //
265 // Trigger jet staff
266 //
267 if(tr.GetNJetThreshold()>0) {
268 // Jet phi/eta
269 Int_t n0 = triggerPosition.GetSize();
270 const TH2F *hpatch = tr.GetJetMatrixE();
271 triggerPosition.Set(n0 + 2);
272 for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);
273 // Add jet ampitudes
274 n0 = triggerAmplitudes.GetSize();
275 triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
276 Double_t *ampJet = tr.GetL1JetThresholds();
277 for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
278 triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
279 }
280 }
0964c2e9 281 esd->AddEMCALTriggerPosition(triggerPosition);
282 esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
85c25c2e 283 // Fill trigger hists
284 AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
aaa3cb7c 285
0e7c6655 286 //########################################
287 //##############Fill CaloCells###############
288 //########################################
aaa3cb7c 289
0e7c6655 290 TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
291 TBranch *branchdig = digitsTree->GetBranch("EMCAL");
292 if (!branchdig) {
0c5b726e 293 AliError("can't get the branch with the EMCAL digits !");
0e7c6655 294 return;
295 }
296 branchdig->SetAddress(&digits);
297 digitsTree->GetEvent(0);
298 Int_t nDigits = digits->GetEntries(), idignew = 0 ;
299 AliDebug(1,Form("%d digits",nDigits));
300
301 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
302 emcCells.CreateContainer(nDigits);
303 emcCells.SetType(AliESDCaloCells::kEMCALCell);
0c5b726e 304 Float_t energy = 0;
0e7c6655 305 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
306 const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
307 if(dig->GetAmp() > 0 ){
0c5b726e 308 energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmp(),dig->GetId());
40164976 309 if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
310 emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());
311 idignew++;
312 }
0e7c6655 313 }
314 }
315 emcCells.SetNumberOfCells(idignew);
316 emcCells.Sort();
317
318 //------------------------------------------------------------
319 //-----------------CLUSTERS-----------------------------
320 //------------------------------------------------------------
321 TObjArray *clusters = new TObjArray(100);
322 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
323 branch->SetAddress(&clusters);
324 clustersTree->GetEvent(0);
325
326 Int_t nClusters = clusters->GetEntries(), nClustersNew=0;
327 AliDebug(1,Form("%d clusters",nClusters));
328 esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters
329
85c25c2e 330
6a0cf740 331 //######################################################
332 //#######################TRACK MATCHING###############
333 //######################################################
334 //Fill list of integers, each one is index of track to which the cluster belongs.
335
336 // step 1 - initialize array of matched track indexes
337 Int_t *matchedTrack = new Int_t[nClusters];
338 for (Int_t iclus = 0; iclus < nClusters; iclus++)
339 matchedTrack[iclus] = -1; // neg. index --> no matched track
340
341 // step 2, change the flag for all matched clusters found in tracks
342 Int_t iemcalMatch = -1;
343 Int_t endtpc = esd->GetNumberOfTracks();
344 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
345 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
346 iemcalMatch = track->GetEMCALcluster();
65f4a419 347 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
6a0cf740 348 }
85c25c2e 349
6a0cf740 350 //########################################
85c25c2e 351 //##############Fill CaloClusters#############
6a0cf740 352 //########################################
35397e76 353 esd->SetNumberOfEMCALClusters(nClusters);
5dee926e 354 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
c47157cd 355 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
8ada0ffe 356 //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
85c60a8e 357 if (Debug()) clust->Print();
a7a5421e 358 // Get information from EMCAL reconstruction points
85c60a8e 359 Float_t xyz[3];
5dee926e 360 TVector3 gpos;
361 clust->GetGlobalPosition(gpos);
35397e76 362 for (Int_t ixyz=0; ixyz<3; ixyz++)
5dee926e 363 xyz[ixyz] = gpos[ixyz];
85c25c2e 364 Float_t elipAxis[2];
365 clust->GetElipsAxis(elipAxis);
35397e76 366 //Create digits lists
367 Int_t cellMult = clust->GetMultiplicity();
368 //TArrayS digiList(digitMult);
369 Float_t *amplFloat = clust->GetEnergiesList();
370 Int_t *digitInts = clust->GetAbsId();
371 TArrayS absIdList(cellMult);
eb972628 372 TArrayD fracList(cellMult);
35397e76 373
374 Int_t newCellMult = 0;
375 for (Int_t iCell=0; iCell<cellMult; iCell++) {
376 if (amplFloat[iCell] > 0) {
377 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
eb972628 378 //Uncomment when unfolding is done
379 //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
380 //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value
381 //else
382 fracList[newCellMult] = 0;
35397e76 383 newCellMult++;
92da3372 384 }
92da3372 385 }
85c25c2e 386
eb972628 387 absIdList.Set(newCellMult);
388 fracList.Set(newCellMult);
389
35397e76 390 if(newCellMult > 0) { // accept cluster if it has some digit
391 nClustersNew++;
65721814 392 //Primaries
7592dfc4 393 Int_t parentMult = 0;
fa42b1f3 394 Int_t *parentList = clust->GetParents(parentMult);
a7a5421e 395 // fills the ESDCaloCluster
35397e76 396 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
397 ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
7592dfc4 398 ec->SetPosition(xyz);
399 ec->SetE(clust->GetEnergy());
40164976 400
401 //Distance to the nearest bad crystal
402 ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower());
403
35397e76 404 ec->SetNCells(newCellMult);
405 //Change type of list from short to ushort
406 UShort_t *newAbsIdList = new UShort_t[newCellMult];
eb972628 407 Double_t *newFracList = new Double_t[newCellMult];
35397e76 408 for(Int_t i = 0; i < newCellMult ; i++) {
409 newAbsIdList[i]=absIdList[i];
eb972628 410 newFracList[i]=fracList[i];
35397e76 411 }
412 ec->SetCellsAbsId(newAbsIdList);
eb972628 413 ec->SetCellsAmplitudeFraction(newFracList);
35397e76 414 ec->SetClusterDisp(clust->GetDispersion());
415 ec->SetClusterChi2(-1); //not yet implemented
416 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
417 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
78902954 418 ec->SetTOF(clust->GetTime()) ; //time-of-fligh
225cd96d 419 ec->SetNExMax(clust->GetNExMax()); //number of local maxima
35397e76 420 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
421 arrayTrackMatched[0]= matchedTrack[iClust];
422 ec->AddTracksMatched(arrayTrackMatched);
423
424 TArrayI arrayParents(parentMult,parentList);
425 ec->AddLabels(arrayParents);
426
6a0cf740 427 // add the cluster to the esd object
eb972628 428 esd->AddCaloCluster(ec);
a7a5421e 429 delete ec;
85d4cbde 430 delete [] newAbsIdList ;
eb972628 431 delete [] newFracList ;
35397e76 432 }
433 } // cycle on clusters
85c25c2e 434
35397e76 435 delete [] matchedTrack;
85c25c2e 436
35397e76 437 esd->SetNumberOfEMCALClusters(nClustersNew);
438 //if(nClustersNew != nClusters)
439 //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
85c25c2e 440
35397e76 441 //Fill ESDCaloCluster with PID weights
85c25c2e 442 AliEMCALPID *pid = new AliEMCALPID;
443 //pid->SetPrintInfo(kTRUE);
444 pid->SetReconstructor(kTRUE);
445 pid->RunPID(esd);
85c25c2e 446 delete pid;
eb972628 447
35397e76 448 delete digits;
85c25c2e 449 delete clusters;
eb972628 450
35397e76 451 // printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew);
0c5b726e 452
453 //Store EMCAL misalignment matrixes
454 FillMisalMatrixes(esd) ;
455
456}
457
458//==================================================================================
459void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
460 //Store EMCAL matrixes in ESD Header
461
462 //Check, if matrixes was already stored
463 for(Int_t sm = 0 ; sm < 12; sm++){
464 if(esd->GetEMCALMatrix(sm)!=0)
465 return ;
466 }
467
468 //Create and store matrixes
469 if(!gGeoManager){
470 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
471 return ;
472 }
473 //Note, that owner of copied marixes will be header
474 char path[255] ;
475 TGeoHMatrix * m ;
476 for(Int_t sm = 0; sm < 12; sm++){
477 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
478 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
479
480 if (gGeoManager->cd(path)){
481 m = gGeoManager->GetCurrentMatrix() ;
482 esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
483 }
484 else{
485 esd->SetEMCALMatrix(NULL,sm) ;
486 }
487 }
488
f6019cda 489}
dc293ae9 490
0c5b726e 491
492
9517d886 493//__________________________________________________________________________
85c25c2e 494void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
495{
496 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
497 if(fgDigitsArr) {
498 // Clear previous digits
499 fgDigitsArr->Delete();
500 delete fgDigitsArr;
501 }
502 // Read the digits from the input tree
503 TBranch *branch = digitsTree->GetBranch("EMCAL");
504 if (!branch) {
505 AliError("can't get the branch with the EMCAL digits !");
506 return;
507 }
508 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
509 branch->SetAddress(&fgDigitsArr);
510 branch->GetEntry(0);
511}
98e9578e 512