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 "AliESDVZERO.h"
57 #include "AliCDBManager.h"
58 #include "AliRunLoader.h"
60 #include "AliEMCALTriggerData.h"
61 #include "AliEMCALTriggerElectronics.h"
62 #include "AliVZEROLoader.h"
64 ClassImp(AliEMCALReconstructor)
66 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
67 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
68 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0; // EMCAL clusterizer class
69 TClonesArray* AliEMCALReconstructor::fgDigitsArr = 0; // shoud read just once at event
70 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
71 //____________________________________________________________________________
72 AliEMCALReconstructor::AliEMCALReconstructor()
73 : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0)
77 fgRawUtils = new AliEMCALRawUtils;
79 //To make sure we match with the geometry in a simulation file,
80 //let's try to get it first. If not, take the default geometry
81 AliRunLoader *rl = AliRunLoader::Instance();
82 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
83 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
85 AliInfo(Form("Using default geometry in reconstruction"));
86 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
89 //Get calibration parameters
92 AliCDBEntry *entry = (AliCDBEntry*)
93 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
94 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
98 AliFatal("Calibration parameters not found in CDB!");
100 //Get calibration parameters
103 AliCDBEntry *entry = (AliCDBEntry*)
104 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
105 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
109 AliFatal("Dead map not found in CDB!");
112 //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
113 fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData);
115 if(!fGeom) AliFatal(Form("Could not get geometry!"));
117 fgTriggerProcessor = new AliEMCALTriggerElectronics();
120 //____________________________________________________________________________
121 AliEMCALReconstructor::~AliEMCALReconstructor()
126 delete fgClusterizer;
127 delete fgTriggerProcessor;
129 AliCodeTimer::Instance()->Print();
132 // //____________________________________________________________________________
133 // void AliEMCALReconstructor::Init()
135 // // Trigger hists - Oct 24, 2007
136 // fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
139 //____________________________________________________________________________
140 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
142 // method called by AliReconstruction;
143 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
144 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
145 // the global tracking.
146 // Works on the current event.
148 AliCodeTimerAuto("",0)
150 ReadDigitsArrayFromTree(digitsTree);
151 fgClusterizer->InitParameters();
152 fgClusterizer->SetOutput(clustersTree);
154 AliEMCALTriggerData* trgData = new AliEMCALTriggerData();
156 Int_t bufferSize = 32000;
158 if (TBranch* triggerBranch = clustersTree->GetBranch("EMTRG"))
159 triggerBranch->SetAddress(&trgData);
161 clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
163 AliVZEROLoader* vzeroLoader = dynamic_cast<AliVZEROLoader*>(AliRunLoader::Instance()->GetDetectorLoader("VZERO"));
169 vzeroLoader->LoadDigits("READ");
170 treeV0 = vzeroLoader->TreeD();
173 TClonesArray *trgDigits = new TClonesArray("AliEMCALRawDigit",1000);
174 TBranch *branchdig = digitsTree->GetBranch("EMTRG");
177 AliError("Can't get the branch with the EMCAL trigger digits !");
181 branchdig->SetAddress(&trgDigits);
182 branchdig->GetEntry(0);
184 fgTriggerProcessor->Digits2Trigger(trgDigits, treeV0, trgData);
188 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
190 fgClusterizer->SetInput(digitsTree);
193 fgClusterizer->Digits2Clusters("deb all") ;
195 fgClusterizer->Digits2Clusters("");
197 fgClusterizer->Clear();
203 vzeroLoader->UnloadDigits();
206 clustersTree->Fill();
211 //____________________________________________________________________________
212 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
215 // Conversion from raw data to
217 // Works on a single-event basis
221 TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
222 TClonesArray *digitsTrg = new TClonesArray("AliEMCALRawDigit", 200);
224 Int_t bufsize = 32000;
225 digitsTree->Branch("EMCAL", &digitsArr, bufsize);
226 digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
228 //must be done here because, in constructor, option is not yet known
229 fgRawUtils->SetOption(GetOption());
231 fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
232 fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
233 fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
234 fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
235 fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
236 fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
237 fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
238 fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
240 fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
251 //____________________________________________________________________________
252 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
253 AliESDEvent* esd) const
255 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
257 // Works on the current event
258 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
261 //######################################################
262 //#########Calculate trigger and set trigger info###########
263 //######################################################
264 // Obsolete, to be changed with new trigger emulator when consensus is achieved about what is stored in ESDs.
266 // tr.SetPatchSize(1); // create 4x4 patches
267 tr.SetSimulation(kFALSE); // Reconstruction mode
268 tr.SetDigitsList(fgDigitsArr);
269 // Get VZERO total multiplicity for jet trigger simulation
270 // The simulation of jey trigger will be incorrect if no VZERO data
272 AliESDVZERO* vZero = esd->GetVZEROData();
274 tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
279 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
280 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
281 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
282 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
284 Int_t iSM2x2 = tr.Get2x2SuperModule();
285 Int_t iSMnxn = tr.GetnxnSuperModule();
286 Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
287 Int_t iModulePhinxn = tr.GetnxnModulePhi();
288 Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
289 Int_t iModuleEtanxn = tr.GetnxnModuleEta();
291 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
292 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
294 TVector3 pos2x2(-1,-1,-1);
295 TVector3 posnxn(-1,-1,-1);
297 Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
298 Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
299 fGeom->GetGlobal(iAbsId2x2, pos2x2);
300 fGeom->GetGlobal(iAbsIdnxn, posnxn);
301 //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
303 TArrayF triggerPosition(6);
304 triggerPosition[0] = pos2x2(0) ;
305 triggerPosition[1] = pos2x2(1) ;
306 triggerPosition[2] = pos2x2(2) ;
307 triggerPosition[3] = posnxn(0) ;
308 triggerPosition[4] = posnxn(1) ;
309 triggerPosition[5] = posnxn(2) ;
310 //printf(" triggerPosition ");
311 //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
313 TArrayF triggerAmplitudes(4);
314 triggerAmplitudes[0] = maxAmp2x2 ;
315 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
316 triggerAmplitudes[2] = maxAmpnxn ;
317 triggerAmplitudes[3] = ampOutOfPatchnxn ;
318 //printf("\n triggerAmplitudes ");
319 //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
325 if(tr.GetNJetThreshold()>0) {
327 Int_t n0 = triggerPosition.GetSize();
328 const TH2F *hpatch = tr.GetJetMatrixE();
329 triggerPosition.Set(n0 + 2);
330 for(Int_t i=0; i<2; i++) triggerPosition[n0+i] = hpatch->GetMean(i+1);
332 n0 = triggerAmplitudes.GetSize();
333 triggerAmplitudes.Set(n0 + tr.GetNJetThreshold());
334 Double_t *ampJet = tr.GetL1JetThresholds();
335 for(Int_t i=0; i<tr.GetNJetThreshold(); i++){
336 triggerAmplitudes[n0 + i] = Float_t(ampJet[i]);
339 esd->AddEMCALTriggerPosition(triggerPosition);
340 esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
341 // Fill trigger hists
342 // AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
344 //########################################
345 //##############Fill CaloCells###############
346 //########################################
348 TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
349 TBranch *branchdig = digitsTree->GetBranch("EMCAL");
351 AliError("can't get the branch with the EMCAL digits !");
354 branchdig->SetAddress(&digits);
355 digitsTree->GetEvent(0);
356 Int_t nDigits = digits->GetEntries(), idignew = 0 ;
357 AliDebug(1,Form("%d digits",nDigits));
359 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
360 emcCells.CreateContainer(nDigits);
361 emcCells.SetType(AliESDCaloCells::kEMCALCell);
363 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
364 const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
365 if(dig->GetAmplitude() > 0 ){
366 energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
367 if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
368 emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());
373 emcCells.SetNumberOfCells(idignew);
376 //------------------------------------------------------------
377 //-----------------CLUSTERS-----------------------------
378 //------------------------------------------------------------
379 TObjArray *clusters = new TObjArray(100);
380 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
381 branch->SetAddress(&clusters);
382 clustersTree->GetEvent(0);
384 Int_t nClusters = clusters->GetEntries(), nClustersNew=0;
385 AliDebug(1,Form("%d clusters",nClusters));
387 //######################################################
388 //#######################TRACK MATCHING###############
389 //######################################################
390 //Fill list of integers, each one is index of track to which the cluster belongs.
392 // step 1 - initialize array of matched track indexes
393 Int_t *matchedTrack = new Int_t[nClusters];
394 for (Int_t iclus = 0; iclus < nClusters; iclus++)
395 matchedTrack[iclus] = -1; // neg. index --> no matched track
397 // step 2, change the flag for all matched clusters found in tracks
398 Int_t iemcalMatch = -1;
399 Int_t endtpc = esd->GetNumberOfTracks();
400 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
401 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
402 iemcalMatch = track->GetEMCALcluster();
403 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
406 //########################################
407 //##############Fill CaloClusters#############
408 //########################################
409 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
410 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
411 //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
412 if (Debug()) clust->Print();
413 // Get information from EMCAL reconstruction points
416 clust->GetGlobalPosition(gpos);
417 for (Int_t ixyz=0; ixyz<3; ixyz++)
418 xyz[ixyz] = gpos[ixyz];
420 clust->GetElipsAxis(elipAxis);
421 //Create digits lists
422 Int_t cellMult = clust->GetMultiplicity();
423 //TArrayS digiList(digitMult);
424 Float_t *amplFloat = clust->GetEnergiesList();
425 Int_t *digitInts = clust->GetAbsId();
426 TArrayS absIdList(cellMult);
427 TArrayD fracList(cellMult);
429 Int_t newCellMult = 0;
430 for (Int_t iCell=0; iCell<cellMult; iCell++) {
431 if (amplFloat[iCell] > 0) {
432 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
433 //Uncomment when unfolding is done
434 //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
435 //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value
437 fracList[newCellMult] = 0;
442 absIdList.Set(newCellMult);
443 fracList.Set(newCellMult);
445 if(newCellMult > 0) { // accept cluster if it has some digit
448 Int_t parentMult = 0;
449 Int_t *parentList = clust->GetParents(parentMult);
450 // fills the ESDCaloCluster
451 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
452 ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
453 ec->SetPosition(xyz);
454 ec->SetE(clust->GetEnergy());
456 //Distance to the nearest bad crystal
457 ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower());
459 ec->SetNCells(newCellMult);
460 //Change type of list from short to ushort
461 UShort_t *newAbsIdList = new UShort_t[newCellMult];
462 Double_t *newFracList = new Double_t[newCellMult];
463 for(Int_t i = 0; i < newCellMult ; i++) {
464 newAbsIdList[i]=absIdList[i];
465 newFracList[i]=fracList[i];
467 ec->SetCellsAbsId(newAbsIdList);
468 ec->SetCellsAmplitudeFraction(newFracList);
469 ec->SetClusterDisp(clust->GetDispersion());
470 ec->SetClusterChi2(-1); //not yet implemented
471 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
472 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
473 ec->SetTOF(clust->GetTime()) ; //time-of-fligh
474 ec->SetNExMax(clust->GetNExMax()); //number of local maxima
475 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
476 arrayTrackMatched[0]= matchedTrack[iClust];
477 ec->AddTracksMatched(arrayTrackMatched);
479 TArrayI arrayParents(parentMult,parentList);
480 ec->AddLabels(arrayParents);
482 // add the cluster to the esd object
483 esd->AddCaloCluster(ec);
485 delete [] newAbsIdList ;
486 delete [] newFracList ;
488 } // cycle on clusters
490 delete [] matchedTrack;
492 //Fill ESDCaloCluster with PID weights
493 AliEMCALPID *pid = new AliEMCALPID;
494 //pid->SetPrintInfo(kTRUE);
495 pid->SetReconstructor(kTRUE);
502 //Store EMCAL misalignment matrixes
503 FillMisalMatrixes(esd) ;
507 //==================================================================================
508 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
509 //Store EMCAL matrixes in ESD Header
511 //Check, if matrixes was already stored
512 for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
513 if(esd->GetEMCALMatrix(sm)!=0)
517 //Create and store matrixes
519 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
522 //Note, that owner of copied marixes will be header
524 TGeoHMatrix * m = 0x0;
525 for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
526 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
527 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
529 if (gGeoManager->CheckPath(path)){
530 gGeoManager->cd(path);
531 m = gGeoManager->GetCurrentMatrix() ;
532 // printf("================================================= \n");
533 // printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
535 esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
536 // printf("================================================= \n");
539 esd->SetEMCALMatrix(NULL,sm) ;
546 //__________________________________________________________________________
547 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
549 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
551 // Clear previous digits
552 fgDigitsArr->Delete();
555 // Read the digits from the input tree
556 TBranch *branch = digitsTree->GetBranch("EMCAL");
558 AliError("can't get the branch with the EMCAL digits !");
561 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
562 branch->SetAddress(&fgDigitsArr);