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.
26 // --- ROOT system ---
27 #include <TClonesArray.h>
28 #include "TGeoManager.h"
29 #include "TGeoMatrix.h"
31 // --- Standard library ---
33 // --- AliRoot header files ---
34 #include "AliEMCALReconstructor.h"
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"
56 #include "AliESDVZERO.h"
57 #include "AliCDBManager.h"
58 #include "AliRunLoader.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"
69 ClassImp(AliEMCALReconstructor)
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)
84 fgRawUtils = new AliEMCALRawUtils;
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();
90 AliEMCAL * emcal = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"));
91 if(emcal) fGeom = emcal->GetGeometry();
95 AliInfo(Form("Using default geometry in reconstruction"));
96 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
99 //Get calibration parameters
102 AliCDBEntry *entry = (AliCDBEntry*)
103 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
104 if (entry) fCalibData = (AliEMCALCalibData*) entry->GetObject();
108 AliFatal("Calibration parameters not found in CDB!");
110 //Get calibration parameters
113 AliCDBEntry *entry = (AliCDBEntry*)
114 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
115 if (entry) fPedestalData = (AliCaloCalibPedestal*) entry->GetObject();
119 AliFatal("Dead map not found in CDB!");
121 if(!fGeom) AliFatal(Form("Could not get geometry!"));
123 AliEMCALTriggerDCSConfigDB* dcsConfigDB = AliEMCALTriggerDCSConfigDB::Instance();
125 const AliEMCALTriggerDCSConfig* dcsConfig = dcsConfigDB->GetTriggerDCSConfig();
127 if (!dcsConfig) AliFatal("No Trigger DCS Configuration from OCDB!");
128 fgTriggerProcessor = new AliEMCALTriggerElectronics( dcsConfig );
130 fTriggerData = new AliEMCALTriggerData();
132 //Init temporary list of digits
133 fgDigitsArr = new TClonesArray("AliEMCALDigit",1000);
134 fgClustersArr = new TObjArray(1000);
135 fgTriggerDigits = new TClonesArray("AliEMCALTriggerRawDigit",1000);
138 //____________________________________________________________________________
139 AliEMCALReconstructor::~AliEMCALReconstructor()
143 if(fGeom) delete fGeom;
145 //No need to delete, recovered from OCDB
146 //if(fCalibData) delete fCalibData;
147 //if(fPedestalData) delete fPedestalData;
150 fgDigitsArr->Clear("C");
155 fgClustersArr->Clear();
156 delete fgClustersArr;
160 fgTriggerDigits->Clear();
161 delete fgTriggerDigits;
164 if(fgRawUtils) delete fgRawUtils;
165 if(fgClusterizer) delete fgClusterizer;
166 if(fgTriggerProcessor) delete fgTriggerProcessor;
168 AliCodeTimer::Instance()->Print();
171 //____________________________________________________________________________
172 void AliEMCALReconstructor::InitClusterizer() const
174 //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
175 Int_t clusterizerType = -1;
176 Int_t eventType = -1;
178 clusterizerType = GetRecParam()->GetClusterizerFlag();
179 eventType = GetRecParam()->GetEventSpecie();
182 AliCDBEntry *entry = (AliCDBEntry*)
183 AliCDBManager::Instance()->Get("EMCAL/Calib/RecoParam");
184 //Get The reco param for the default event specie
186 AliEMCALRecParam *recParam = (AliEMCALRecParam*)((TObjArray *) entry->GetObject())->At(0);
187 if(recParam) clusterizerType = recParam->GetClusterizerFlag();
191 //Check if clusterizer previously set corresponds to what is needed for this event type
193 if(eventType!=AliRecoParam::kCalib){
194 //printf("ReCreate clusterizer? Clusterizer set <%d>, Clusterizer in use <%s>\n",
195 // clusterizerType, fgClusterizer->Version());
197 if (clusterizerType == AliEMCALRecParam::kClusterizerv1 && !strcmp(fgClusterizer->Version(),"clu-v1")) return;
199 else if(clusterizerType == AliEMCALRecParam::kClusterizerNxN && !strcmp(fgClusterizer->Version(),"clu-NxN")) return;
201 //Need to create new clusterizer, the one set previously is not the correct one
202 delete fgClusterizer;
207 if (clusterizerType == AliEMCALRecParam::kClusterizerv1)
209 fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData);
213 fgClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData,fPedestalData);
217 //____________________________________________________________________________
218 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
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.
226 AliCodeTimerAuto("",0)
228 //Get input digits and put them in fgDigitsArr, clear the list before
229 ReadDigitsArrayFromTree(digitsTree);
233 fgClusterizer->InitParameters();
234 fgClusterizer->SetOutput(clustersTree);
236 //Skip clusterization of LED events
237 if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
239 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
241 fgClusterizer->SetInput(digitsTree);
243 //fgClusterizer->Digits2Clusters("deb all") ; //For debugging
244 fgClusterizer->Digits2Clusters("");
246 fgClusterizer->Clear();
248 }//digits array exists and has somethind
251 clustersTree->Fill();
254 //____________________________________________________________________________
255 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
258 // Conversion from raw data to
260 // Works on a single-event basis
264 fTriggerData->SetMode(1);
266 if(fgDigitsArr) fgDigitsArr->Clear("C");
268 TClonesArray *digitsTrg = new TClonesArray("AliEMCALTriggerRawDigit", 32 * 96);
270 Int_t bufsize = 32000;
271 digitsTree->Branch("EMCAL", &fgDigitsArr, bufsize);
272 digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
274 //Skip calibration events do the rest
275 Bool_t doFit = kTRUE;
276 if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
278 //must be done here because, in constructor, option is not yet known
279 fgRawUtils->SetOption(GetOption());
281 fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
282 fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
283 fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
284 fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
285 fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
286 fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
287 fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
288 fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
289 fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
290 fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
292 fgRawUtils->Raw2Digits(rawReader,fgDigitsArr,fPedestalData,digitsTrg,fTriggerData);
293 }//skip calibration event
295 AliDebug(1," Calibration Event, skip!");
305 //____________________________________________________________________________
306 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
307 AliESDEvent* esd) const
309 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
311 // Works on the current event
312 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
315 //########################################
317 //########################################
319 Int_t v0M[2] = {0, 0};
321 AliESDVZERO* esdV0 = esd->GetVZEROData();
325 for (Int_t i = 0; i < 32; i++)
327 v0M[0] += (Int_t)esdV0->GetAdcV0C(i);
328 v0M[1] += (Int_t)esdV0->GetAdcV0A(i);
333 AliWarning("Cannot retrieve V0 ESD! Run w/ null V0 charges");
336 if (fgTriggerDigits) fgTriggerDigits->Clear();
338 TBranch *branchtrg = digitsTree->GetBranch("EMTRG");
342 AliError("Can't get the branch with the EMCAL trigger digits!");
346 branchtrg->SetAddress(&fgTriggerDigits);
347 branchtrg->GetEntry(0);
349 // Note: fgTriggerProcessor reset done at the end of this method
350 fgTriggerProcessor->Digits2Trigger(fgTriggerDigits, v0M, fTriggerData);
353 AliESDCaloTrigger* trgESD = esd->GetCaloTrigger("EMCAL");
357 trgESD->Allocate(fgTriggerDigits->GetEntriesFast());
359 for (Int_t i = 0; i < fgTriggerDigits->GetEntriesFast(); i++)
361 AliEMCALTriggerRawDigit* rdig = (AliEMCALTriggerRawDigit*)fgTriggerDigits->At(i);
364 if (fGeom->GetPositionInEMCALFromAbsFastORIndex(rdig->GetId(), px, py))
366 Int_t a = -1, t = -1, times[10];
368 rdig->GetMaximum(a, t);
369 rdig->GetL0Times(times);
371 trgESD->Add(px, py, a, t, times, rdig->GetNL0Times(), rdig->GetL1TimeSum(), rdig->GetTriggerBits());
375 trgESD->SetL1Threshold(0, fTriggerData->GetL1GammaThreshold());
377 trgESD->SetL1Threshold(1, fTriggerData->GetL1JetThreshold() );
381 fTriggerData->Reset();
383 //########################################
384 //##############Fill CaloCells###############
385 //########################################
387 //Get input digits and put them in fgDigitsArr, clear the list before
388 ReadDigitsArrayFromTree(digitsTree);
390 Int_t nDigits = fgDigitsArr->GetEntries(), idignew = 0 ;
391 AliDebug(1,Form("%d digits",nDigits));
393 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
394 emcCells.CreateContainer(nDigits);
395 emcCells.SetType(AliVCaloCells::kEMCALCell);
397 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
398 const AliEMCALDigit * dig = (const AliEMCALDigit*)fgDigitsArr->At(idig);
399 if(dig->GetAmplitude() > 0 ){
400 energy = fgClusterizer->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
401 if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them
402 emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());
407 emcCells.SetNumberOfCells(idignew);
410 //------------------------------------------------------------
411 //-----------------CLUSTERS-----------------------------
412 //------------------------------------------------------------
413 clustersTree->SetBranchStatus("*",0); //disable all branches
414 clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
415 if(fgClustersArr) fgClustersArr->Clear();
416 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
417 branch->SetAddress(&fgClustersArr);
419 //clustersTree->GetEvent(0);
421 Int_t nClusters = fgClustersArr->GetEntries(), nClustersNew=0;
422 AliDebug(1,Form("%d clusters",nClusters));
424 //######################################################
425 //#######################TRACK MATCHING###############
426 //######################################################
427 //Fill list of integers, each one is index of track to which the cluster belongs.
429 // step 1 - initialize array of matched track indexes
430 Int_t *matchedTrack = new Int_t[nClusters];
431 for (Int_t iclus = 0; iclus < nClusters; iclus++)
432 matchedTrack[iclus] = -1; // neg. index --> no matched track
434 // step 2, change the flag for all matched clusters found in tracks
435 Int_t iemcalMatch = -1;
436 Int_t endtpc = esd->GetNumberOfTracks();
437 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
438 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
439 iemcalMatch = track->GetEMCALcluster();
440 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
443 //########################################
444 //##############Fill CaloClusters#############
445 //########################################
446 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
447 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)fgClustersArr->At(iClust);
448 //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
449 // clust->Print(); //For debugging
450 // Get information from EMCAL reconstruction points
453 clust->GetGlobalPosition(gpos);
454 for (Int_t ixyz=0; ixyz<3; ixyz++)
455 xyz[ixyz] = gpos[ixyz];
457 clust->GetElipsAxis(elipAxis);
458 //Create digits lists
459 Int_t cellMult = clust->GetMultiplicity();
460 //TArrayS digiList(digitMult);
461 Float_t *amplFloat = clust->GetEnergiesList();
462 Int_t *digitInts = clust->GetAbsId();
463 TArrayS absIdList(cellMult);
464 TArrayD fracList(cellMult);
466 Int_t newCellMult = 0;
467 for (Int_t iCell=0; iCell<cellMult; iCell++) {
468 if (amplFloat[iCell] > 0) {
469 absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
471 if(emcCells.GetCellAmplitude(digitInts[iCell])>0 && GetRecParam()->GetUnfold()){
472 fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell]));//get cell calibration value
476 fracList[newCellMult] = 0;
482 absIdList.Set(newCellMult);
483 fracList.Set(newCellMult);
485 if(newCellMult > 0) { // accept cluster if it has some digit
488 Int_t parentMult = 0;
489 Int_t *parentList = clust->GetParents(parentMult);
490 // fills the ESDCaloCluster
491 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
492 ec->SetType(AliVCluster::kEMCALClusterv1);
493 ec->SetPosition(xyz);
494 ec->SetE(clust->GetEnergy());
496 //Distance to the nearest bad crystal
497 ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower());
499 ec->SetNCells(newCellMult);
500 //Change type of list from short to ushort
501 UShort_t *newAbsIdList = new UShort_t[newCellMult];
502 Double_t *newFracList = new Double_t[newCellMult];
503 for(Int_t i = 0; i < newCellMult ; i++) {
504 newAbsIdList[i]=absIdList[i];
505 newFracList[i] =fracList[i];
507 ec->SetCellsAbsId(newAbsIdList);
508 ec->SetCellsAmplitudeFraction(newFracList);
509 ec->SetDispersion(clust->GetDispersion());
510 ec->SetChi2(-1); //not yet implemented
511 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
512 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
513 ec->SetTOF(clust->GetTime()) ; //time-of-fligh
514 ec->SetNExMax(clust->GetNExMax()); //number of local maxima
515 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
516 arrayTrackMatched[0]= matchedTrack[iClust];
517 ec->AddTracksMatched(arrayTrackMatched);
519 TArrayI arrayParents(parentMult,parentList);
520 ec->AddLabels(arrayParents);
522 // add the cluster to the esd object
523 esd->AddCaloCluster(ec);
525 delete [] newAbsIdList ;
526 delete [] newFracList ;
528 } // cycle on clusters
530 delete [] matchedTrack;
532 //Fill ESDCaloCluster with PID weights
533 AliEMCALPID *pid = new AliEMCALPID;
534 //pid->SetPrintInfo(kTRUE);
535 pid->SetReconstructor(kTRUE);
539 //Store EMCAL misalignment matrixes
540 FillMisalMatrixes(esd) ;
544 //==================================================================================
545 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
546 //Store EMCAL matrixes in ESD Header
548 //Check, if matrixes was already stored
549 for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
550 if(esd->GetEMCALMatrix(sm)!=0)
554 //Create and store matrixes
556 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
559 //Note, that owner of copied marixes will be header
560 const Int_t bufsize = 255;
562 TGeoHMatrix * m = 0x0;
563 for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
564 snprintf(path,bufsize,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
565 if(sm >= 10) snprintf(path,bufsize,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
567 if (gGeoManager->CheckPath(path)){
568 gGeoManager->cd(path);
569 m = gGeoManager->GetCurrentMatrix() ;
570 // printf("================================================= \n");
571 // printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
573 esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
574 // printf("================================================= \n");
577 esd->SetEMCALMatrix(NULL,sm) ;
582 //__________________________________________________________________________
583 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
585 // Read the digits from the input tree
586 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
588 // Clear previous digits in the list
590 fgDigitsArr->Clear("C");
593 // It should not happen, but just in case ...
594 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
597 // Read the digits from the input tree
598 TBranch *branch = digitsTree->GetBranch("EMCAL");
600 AliError("can't get the branch with the EMCAL digits !");
604 branch->SetAddress(&fgDigitsArr);