1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
18 //_________________________________________________________________________
20 //-- Yves Schutz (SUBATECH)
21 // Reconstruction class. Redesigned from the old AliReconstructionner class and
22 // derived from STEER/AliReconstructor.
24 //-- Aleksei Pavlinov : added staf for EMCAL jet trigger 9Apr 25, 2008)
25 // : fgDigitsArr should read just once at event
27 // --- ROOT system ---
29 #include <TClonesArray.h>
31 #include "TGeoManager.h"
32 #include "TGeoMatrix.h"
34 // --- Standard library ---
36 // --- AliRoot header files ---
37 #include "AliEMCALReconstructor.h"
39 #include "AliCodeTimer.h"
40 #include "AliESDEvent.h"
41 #include "AliESDCaloCluster.h"
42 #include "AliESDCaloCells.h"
43 #include "AliESDtrack.h"
44 #include "AliEMCALLoader.h"
45 #include "AliEMCALRawUtils.h"
46 #include "AliEMCALDigit.h"
47 #include "AliEMCALClusterizerv1.h"
48 #include "AliEMCALRecPoint.h"
49 #include "AliEMCALPID.h"
50 #include "AliEMCALTrigger.h"
51 #include "AliRawReader.h"
52 #include "AliCDBEntry.h"
53 #include "AliCDBManager.h"
54 #include "AliEMCALGeometry.h"
56 #include "AliEMCALHistoUtilities.h"
57 #include "AliESDVZERO.h"
58 #include "AliCDBManager.h"
59 #include "AliRunLoader.h"
62 ClassImp(AliEMCALReconstructor)
64 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
65 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
66 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0; // EMCAL clusterizer class
67 TClonesArray* AliEMCALReconstructor::fgDigitsArr = 0; // shoud read just once at event
68 //____________________________________________________________________________
69 AliEMCALReconstructor::AliEMCALReconstructor()
70 : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0)
74 fgRawUtils = new AliEMCALRawUtils;
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
78 AliRunLoader *rl = AliRunLoader::Instance();
79 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
80 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
82 AliInfo(Form("Using default geometry in reconstruction"));
83 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
86 //Get calibration parameters
89 AliCDBEntry *entry = (AliCDBEntry*)
90 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
91 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
95 AliFatal("Calibration parameters not found in CDB!");
97 //Get calibration parameters
100 AliCDBEntry *entry = (AliCDBEntry*)
101 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
102 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
106 AliFatal("Dead map not found in CDB!");
109 //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
110 fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData);
112 if(!fGeom) AliFatal(Form("Could not get geometry!"));
116 //____________________________________________________________________________
117 AliEMCALReconstructor::~AliEMCALReconstructor()
122 delete fgClusterizer;
124 AliCodeTimer::Instance()->Print();
127 //____________________________________________________________________________
128 void AliEMCALReconstructor::Init()
130 // Trigger hists - Oct 24, 2007
131 fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
134 //____________________________________________________________________________
135 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
137 // method called by AliReconstruction;
138 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
139 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
140 // the global tracking.
141 // Works on the current event.
143 AliCodeTimerAuto("",0)
145 ReadDigitsArrayFromTree(digitsTree);
146 fgClusterizer->InitParameters();
147 fgClusterizer->SetOutput(clustersTree);
149 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
151 fgClusterizer->SetInput(digitsTree);
154 fgClusterizer->Digits2Clusters("deb all") ;
156 fgClusterizer->Digits2Clusters("");
158 fgClusterizer->Clear();
164 //____________________________________________________________________________
165 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
168 // Conversion from raw data to
170 // Works on a single-event basis
174 TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
175 Int_t bufsize = 32000;
176 digitsTree->Branch("EMCAL", &digitsArr, bufsize);
178 //must be done here because, in constructor, option is not yet known
179 fgRawUtils->SetOption(GetOption());
181 fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
182 fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
183 fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
184 fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
185 fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
186 fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
187 fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
189 fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData);
198 //____________________________________________________________________________
199 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
200 AliESDEvent* esd) const
202 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
204 // Works on the current event
205 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
208 //######################################################
209 //#########Calculate trigger and set trigger info###########
210 //######################################################
213 // tr.SetPatchSize(1); // create 4x4 patches
214 tr.SetSimulation(kFALSE); // Reconstruction mode
215 tr.SetDigitsList(fgDigitsArr);
216 // Get VZERO total multiplicity for jet trigger simulation
217 // The simulation of jey trigger will be incorrect if no VZERO data
219 AliESDVZERO* vZero = esd->GetVZEROData();
221 tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
226 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
227 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
228 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
229 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
231 Int_t iSM2x2 = tr.Get2x2SuperModule();
232 Int_t iSMnxn = tr.GetnxnSuperModule();
233 Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
234 Int_t iModulePhinxn = tr.GetnxnModulePhi();
235 Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
236 Int_t iModuleEtanxn = tr.GetnxnModuleEta();
238 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
239 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
241 TVector3 pos2x2(-1,-1,-1);
242 TVector3 posnxn(-1,-1,-1);
244 Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
245 Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
246 fGeom->GetGlobal(iAbsId2x2, pos2x2);
247 fGeom->GetGlobal(iAbsIdnxn, posnxn);
248 //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
250 TArrayF triggerPosition(6);
251 triggerPosition[0] = pos2x2(0) ;
252 triggerPosition[1] = pos2x2(1) ;
253 triggerPosition[2] = pos2x2(2) ;
254 triggerPosition[3] = posnxn(0) ;
255 triggerPosition[4] = posnxn(1) ;
256 triggerPosition[5] = posnxn(2) ;
257 //printf(" triggerPosition ");
258 //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
260 TArrayF triggerAmplitudes(4);
261 triggerAmplitudes[0] = maxAmp2x2 ;
262 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
263 triggerAmplitudes[2] = maxAmpnxn ;
264 triggerAmplitudes[3] = ampOutOfPatchnxn ;
265 //printf("\n triggerAmplitudes ");
266 //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
272 if(tr.GetNJetThreshold()>0) {
274 Int_t n0 = triggerPosition.GetSize();
275 const TH2F *hpatch = tr.GetJetMatrixE();
276 triggerPosition.Set(n0 + 2);
277 for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);
279 n0 = triggerAmplitudes.GetSize();
280 triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
281 Double_t *ampJet = tr.GetL1JetThresholds();
282 for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
283 triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
286 esd->AddEMCALTriggerPosition(triggerPosition);
287 esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
288 // Fill trigger hists
289 AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
291 //########################################
292 //##############Fill CaloCells###############
293 //########################################
295 TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
296 TBranch *branchdig = digitsTree->GetBranch("EMCAL");
298 AliError("can't get the branch with the EMCAL digits !");
301 branchdig->SetAddress(&digits);
302 digitsTree->GetEvent(0);
303 Int_t nDigits = digits->GetEntries(), idignew = 0 ;
304 AliDebug(1,Form("%d digits",nDigits));
306 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
307 emcCells.CreateContainer(nDigits);
308 emcCells.SetType(AliESDCaloCells::kEMCALCell);
310 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
311 const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
312 if(dig->GetAmp() > 0 ){
313 energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmp(),dig->GetId());
314 if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
315 emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());
320 emcCells.SetNumberOfCells(idignew);
323 //------------------------------------------------------------
324 //-----------------CLUSTERS-----------------------------
325 //------------------------------------------------------------
326 TObjArray *clusters = new TObjArray(100);
327 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
328 branch->SetAddress(&clusters);
329 clustersTree->GetEvent(0);
331 Int_t nClusters = clusters->GetEntries(), nClustersNew=0;
332 AliDebug(1,Form("%d clusters",nClusters));
333 esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters
336 //######################################################
337 //#######################TRACK MATCHING###############
338 //######################################################
339 //Fill list of integers, each one is index of track to which the cluster belongs.
341 // step 1 - initialize array of matched track indexes
342 Int_t *matchedTrack = new Int_t[nClusters];
343 for (Int_t iclus = 0; iclus < nClusters; iclus++)
344 matchedTrack[iclus] = -1; // neg. index --> no matched track
346 // step 2, change the flag for all matched clusters found in tracks
347 Int_t iemcalMatch = -1;
348 Int_t endtpc = esd->GetNumberOfTracks();
349 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
350 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
351 iemcalMatch = track->GetEMCALcluster();
352 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
355 //########################################
356 //##############Fill CaloClusters#############
357 //########################################
358 esd->SetNumberOfEMCALClusters(nClusters);
359 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
360 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
361 //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
362 if (Debug()) clust->Print();
363 // Get information from EMCAL reconstruction points
366 clust->GetGlobalPosition(gpos);
367 for (Int_t ixyz=0; ixyz<3; ixyz++)
368 xyz[ixyz] = gpos[ixyz];
370 clust->GetElipsAxis(elipAxis);
371 //Create digits lists
372 Int_t cellMult = clust->GetMultiplicity();
373 //TArrayS digiList(digitMult);
374 Float_t *amplFloat = clust->GetEnergiesList();
375 Int_t *digitInts = clust->GetAbsId();
376 TArrayS absIdList(cellMult);
377 TArrayD fracList(cellMult);
379 Int_t newCellMult = 0;
380 for (Int_t iCell=0; iCell<cellMult; iCell++) {
381 if (amplFloat[iCell] > 0) {
382 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
383 //Uncomment when unfolding is done
384 //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
385 //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value
387 fracList[newCellMult] = 0;
392 absIdList.Set(newCellMult);
393 fracList.Set(newCellMult);
395 if(newCellMult > 0) { // accept cluster if it has some digit
398 Int_t parentMult = 0;
399 Int_t *parentList = clust->GetParents(parentMult);
400 // fills the ESDCaloCluster
401 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
402 ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
403 ec->SetPosition(xyz);
404 ec->SetE(clust->GetEnergy());
406 //Distance to the nearest bad crystal
407 ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower());
409 ec->SetNCells(newCellMult);
410 //Change type of list from short to ushort
411 UShort_t *newAbsIdList = new UShort_t[newCellMult];
412 Double_t *newFracList = new Double_t[newCellMult];
413 for(Int_t i = 0; i < newCellMult ; i++) {
414 newAbsIdList[i]=absIdList[i];
415 newFracList[i]=fracList[i];
417 ec->SetCellsAbsId(newAbsIdList);
418 ec->SetCellsAmplitudeFraction(newFracList);
419 ec->SetClusterDisp(clust->GetDispersion());
420 ec->SetClusterChi2(-1); //not yet implemented
421 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
422 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
423 ec->SetTOF(clust->GetTime()) ; //time-of-fligh
424 ec->SetNExMax(clust->GetNExMax()); //number of local maxima
425 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
426 arrayTrackMatched[0]= matchedTrack[iClust];
427 ec->AddTracksMatched(arrayTrackMatched);
429 TArrayI arrayParents(parentMult,parentList);
430 ec->AddLabels(arrayParents);
432 // add the cluster to the esd object
433 esd->AddCaloCluster(ec);
435 delete [] newAbsIdList ;
436 delete [] newFracList ;
438 } // cycle on clusters
440 delete [] matchedTrack;
442 esd->SetNumberOfEMCALClusters(nClustersNew);
443 //if(nClustersNew != nClusters)
444 //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
446 //Fill ESDCaloCluster with PID weights
447 AliEMCALPID *pid = new AliEMCALPID;
448 //pid->SetPrintInfo(kTRUE);
449 pid->SetReconstructor(kTRUE);
456 // printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew);
458 //Store EMCAL misalignment matrixes
459 FillMisalMatrixes(esd) ;
463 //==================================================================================
464 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
465 //Store EMCAL matrixes in ESD Header
467 //Check, if matrixes was already stored
468 for(Int_t sm = 0 ; sm < 12; sm++){
469 if(esd->GetEMCALMatrix(sm)!=0)
473 //Create and store matrixes
475 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
478 //Note, that owner of copied marixes will be header
480 TGeoHMatrix * m = 0x0;
481 for(Int_t sm = 0; sm < 12; sm++){
482 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
483 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
485 if (gGeoManager->CheckPath(path)){
486 m = gGeoManager->GetCurrentMatrix() ;
487 esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
490 esd->SetEMCALMatrix(NULL,sm) ;
497 //__________________________________________________________________________
498 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
500 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
502 // Clear previous digits
503 fgDigitsArr->Delete();
506 // Read the digits from the input tree
507 TBranch *branch = digitsTree->GetBranch("EMCAL");
509 AliError("can't get the branch with the EMCAL digits !");
512 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
513 branch->SetAddress(&fgDigitsArr);