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 // --- Standard library ---
33 // --- AliRoot header files ---
34 #include "AliEMCALReconstructor.h"
36 #include "AliCodeTimer.h"
37 #include "AliESDEvent.h"
38 #include "AliESDCaloCluster.h"
39 #include "AliESDCaloCells.h"
40 #include "AliESDtrack.h"
41 #include "AliEMCALLoader.h"
42 #include "AliEMCALRawUtils.h"
43 #include "AliEMCALDigit.h"
44 #include "AliEMCALClusterizerv1.h"
45 #include "AliEMCALRecPoint.h"
46 #include "AliEMCALPID.h"
47 #include "AliEMCALTrigger.h"
48 #include "AliRawReader.h"
49 #include "AliCDBEntry.h"
50 #include "AliCDBManager.h"
51 #include "AliEMCALRecParam.h"
52 #include "AliEMCALGeometry.h"
54 #include "AliEMCALHistoUtilities.h"
55 #include "AliESDVZERO.h"
57 #include "AliRunLoader.h"
60 ClassImp(AliEMCALReconstructor)
62 AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0; // EMCAL rec. parameters
63 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0; // EMCAL raw utilities class
64 TClonesArray* AliEMCALReconstructor::fgDigitsArr = 0; // shoud read just once at event
66 //____________________________________________________________________________
67 AliEMCALReconstructor::AliEMCALReconstructor()
68 : fDebug(kFALSE), fList(0), fGeom(0)
73 fgRawUtils = new AliEMCALRawUtils;
75 //To make sure we match with the geometry in a simulation file,
76 //let's try to get it first. If not, take the default geometry
77 AliRunLoader *rl = AliRunLoader::GetRunLoader();
78 if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
79 fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
81 AliInfo(Form("Using default geometry in reconstruction"));
82 fGeom = AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
85 if(!fGeom) AliFatal(Form("Could not get geometry!"));
89 //____________________________________________________________________________
90 AliEMCALReconstructor::AliEMCALReconstructor(const AliEMCALReconstructor & rec)
91 : AliReconstructor(rec),
99 //____________________________________________________________________________
100 AliEMCALReconstructor::~AliEMCALReconstructor()
104 AliCodeTimer::Instance()->Print();
107 //____________________________________________________________________________
108 void AliEMCALReconstructor::Init()
110 // Trigger hists - Oct 24, 2007
111 fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
114 //____________________________________________________________________________
115 void AliEMCALReconstructor::InitRecParam() const
117 // Check if the instance of AliEMCALRecParam exists,
118 // if not, get it from OCDB if available, otherwise create a default one
120 if (!fgkRecParam && (AliCDBManager::Instance()->IsDefaultStorageSet())) {
121 AliCDBEntry *entry = (AliCDBEntry*)
122 AliCDBManager::Instance()->Get("EMCAL/Config/RecParam");
123 if (entry) fgkRecParam = (AliEMCALRecParam*) entry->GetObject();
127 AliWarning("The Reconstruction parameters for EMCAL nonitialized - Used default one");
128 fgkRecParam = new AliEMCALRecParam;
132 //____________________________________________________________________________
133 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
135 // method called by AliReconstruction;
136 // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
137 // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by
138 // the global tracking.
139 // Works on the current event.
143 ReadDigitsArrayFromTree(digitsTree);
145 if(fgDigitsArr && fgDigitsArr->GetEntries()) {
147 AliEMCALClusterizerv1 clu;
148 clu.SetInput(digitsTree);
149 clu.SetOutput(clustersTree);
151 clu.Digits2Clusters("deb all") ;
153 clu.Digits2Clusters("") ;
158 //____________________________________________________________________________
159 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
162 // Conversion from raw data to
164 // Works on a single-event basis
168 TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",100);
169 Int_t bufsize = 32000;
170 digitsTree->Branch("EMCAL", &digitsArr, bufsize);
172 //must be done here because, in constructor, option is not yet known
173 fgRawUtils->SetOption(GetOption());
175 fgRawUtils->SetRawFormatHighLowGainFactor(fgkRecParam->GetHighLowGainFactor());
176 fgRawUtils->SetRawFormatOrder(fgkRecParam->GetOrderParameter());
177 fgRawUtils->SetRawFormatTau(fgkRecParam->GetTau());
178 fgRawUtils->SetNoiseThreshold(fgkRecParam->GetNoiseThreshold());
179 fgRawUtils->SetNPedSamples(fgkRecParam->GetNPedSamples());
181 fgRawUtils->Raw2Digits(rawReader,digitsArr);
190 //____________________________________________________________________________
191 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree,
192 AliESDEvent* esd) const
194 // Called by AliReconstruct after Reconstruct() and global tracking and vertexing
196 // Works on the current event
197 // printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
199 const double timeScale = 1.e+11; // transition constant from sec to 0.01 ns
201 //######################################################
202 //#########Calculate trigger and set trigger info###########
203 //######################################################
206 // tr.SetPatchSize(1); // create 4x4 patches
207 tr.SetSimulation(kFALSE); // Reconstruction mode
208 tr.SetDigitsList(fgDigitsArr);
209 // Get VZERO total multiplicity for jet trigger simulation
210 // The simulation of jey trigger will be incorrect if no VZERO data
212 AliESDVZERO* vZero = esd->GetVZEROData();
214 tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
219 Float_t maxAmp2x2 = tr.Get2x2MaxAmplitude();
220 Float_t maxAmpnxn = tr.GetnxnMaxAmplitude();
221 Float_t ampOutOfPatch2x2 = tr.Get2x2AmpOutOfPatch() ;
222 Float_t ampOutOfPatchnxn = tr.GetnxnAmpOutOfPatch() ;
224 Int_t iSM2x2 = tr.Get2x2SuperModule();
225 Int_t iSMnxn = tr.GetnxnSuperModule();
226 Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
227 Int_t iModulePhinxn = tr.GetnxnModulePhi();
228 Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
229 Int_t iModuleEtanxn = tr.GetnxnModuleEta();
231 AliDebug(2, Form("Trigger 2x2 max amp %f, out amp %f, SM %d, iphi %d ieta %d", maxAmp2x2, ampOutOfPatch2x2, iSM2x2,iModulePhi2x2, iModuleEta2x2));
232 AliDebug(2, Form("Trigger 4x4 max amp %f , out amp %f, SM %d, iphi %d, ieta %d", maxAmpnxn, ampOutOfPatchnxn, iSMnxn,iModulePhinxn, iModuleEtanxn));
234 TVector3 pos2x2(-1,-1,-1);
235 TVector3 posnxn(-1,-1,-1);
237 Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
238 Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
239 fGeom->GetGlobal(iAbsId2x2, pos2x2);
240 fGeom->GetGlobal(iAbsIdnxn, posnxn);
241 //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
243 TArrayF triggerPosition(6);
244 triggerPosition[0] = pos2x2(0) ;
245 triggerPosition[1] = pos2x2(1) ;
246 triggerPosition[2] = pos2x2(2) ;
247 triggerPosition[3] = posnxn(0) ;
248 triggerPosition[4] = posnxn(1) ;
249 triggerPosition[5] = posnxn(2) ;
250 //printf(" triggerPosition ");
251 //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
253 TArrayF triggerAmplitudes(4);
254 triggerAmplitudes[0] = maxAmp2x2 ;
255 triggerAmplitudes[1] = ampOutOfPatch2x2 ;
256 triggerAmplitudes[2] = maxAmpnxn ;
257 triggerAmplitudes[3] = ampOutOfPatchnxn ;
258 //printf("\n triggerAmplitudes ");
259 //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
263 esd->AddEMCALTriggerPosition(triggerPosition);
264 esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
265 // Fill trigger hists
266 AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
268 //########################################
269 //##############Fill CaloCells###############
270 //########################################
272 TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
273 TBranch *branchdig = digitsTree->GetBranch("EMCAL");
275 AliError("can't get the branch with the PHOS digits !");
278 branchdig->SetAddress(&digits);
279 digitsTree->GetEvent(0);
280 Int_t nDigits = digits->GetEntries(), idignew = 0 ;
281 AliDebug(1,Form("%d digits",nDigits));
283 AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
284 emcCells.CreateContainer(nDigits);
285 emcCells.SetType(AliESDCaloCells::kEMCALCell);
286 for (Int_t idig = 0 ; idig < nDigits ; idig++) {
287 const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
288 if(dig->GetAmp() > 0 ){
289 emcCells.SetCell(idignew,dig->GetId(),dig->GetAmp(), dig->GetTime());
293 emcCells.SetNumberOfCells(idignew);
296 //------------------------------------------------------------
297 //-----------------CLUSTERS-----------------------------
298 //------------------------------------------------------------
299 TObjArray *clusters = new TObjArray(100);
300 TBranch *branch = clustersTree->GetBranch("EMCALECARP");
301 branch->SetAddress(&clusters);
302 clustersTree->GetEvent(0);
304 Int_t nClusters = clusters->GetEntries(), nClustersNew=0;
305 AliDebug(1,Form("%d clusters",nClusters));
306 esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters
309 //######################################################
310 //#######################TRACK MATCHING###############
311 //######################################################
312 //Fill list of integers, each one is index of track to which the cluster belongs.
314 // step 1 - initialize array of matched track indexes
315 Int_t *matchedTrack = new Int_t[nClusters];
316 for (Int_t iclus = 0; iclus < nClusters; iclus++)
317 matchedTrack[iclus] = -1; // neg. index --> no matched track
319 // step 2, change the flag for all matched clusters found in tracks
320 Int_t iemcalMatch = -1;
321 Int_t endtpc = esd->GetNumberOfTracks();
322 for (Int_t itrack = 0; itrack < endtpc; itrack++) {
323 AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
324 iemcalMatch = track->GetEMCALcluster();
325 if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
328 //########################################
329 //##############Fill CaloClusters#############
330 //########################################
332 for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
333 const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
334 //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
335 if (Debug()) clust->Print();
336 // Get information from EMCAL reconstruction points
339 clust->GetGlobalPosition(gpos);
340 for (Int_t ixyz=0; ixyz<3; ixyz++)
341 xyz[ixyz] = gpos[ixyz];
343 Int_t digitMult = clust->GetMultiplicity();
344 TArrayS amplList(digitMult);
345 TArrayS timeList(digitMult);
346 TArrayS digiList(digitMult);
347 Float_t *amplFloat = clust->GetEnergiesList();
348 Float_t *timeFloat = clust->GetTimeList();
349 Int_t *digitInts = clust->GetAbsId();
351 clust->GetElipsAxis(elipAxis);
353 // Convert Float_t* and Int_t* to Short_t* to save memory
354 // Problem : we should recalculate a cluster characteristics when discard digit(s)
355 Int_t newdigitMult = 0;
356 for (Int_t iDigit=0; iDigit<digitMult; iDigit++) {
357 if (amplFloat[iDigit] > 0) {
358 amplList[newdigitMult] = (UShort_t)(amplFloat[iDigit]*500);
359 // Time in units of 0.01 ns = 10 ps
360 if(timeFloat[iDigit] < 65536./timeScale)
361 timeList[newdigitMult] = (UShort_t)(timeFloat[iDigit]*timeScale);
363 timeList[newdigitMult] = 65535;
364 digiList[newdigitMult] = (UShort_t)(digitInts[iDigit]);
367 else if (clust->GetClusterType() != AliESDCaloCluster::kEMCALPseudoCluster)
368 Warning("FillESD()","Negative or 0 digit amplitude in cluster");
371 if(newdigitMult > 0) { // accept cluster if it has some digit
374 amplList.Set(newdigitMult);
375 timeList.Set(newdigitMult);
376 digiList.Set(newdigitMult);
379 Int_t parentMult = 0;
380 Int_t *parentList = clust->GetParents(parentMult);
382 // fills the ESDCaloCluster
383 AliESDCaloCluster * ec = new AliESDCaloCluster() ;
384 ec->SetClusterType(clust->GetClusterType());
385 ec->SetPosition(xyz);
386 ec->SetE(clust->GetEnergy());
387 ec->AddDigitAmplitude(amplList);
388 ec->AddDigitTime(timeList);
389 ec->AddDigitIndex(digiList);
391 if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1){
393 ec->SetClusterDisp(clust->GetDispersion());
394 ec->SetClusterChi2(-1); //not yet implemented
395 ec->SetM02(elipAxis[0]*elipAxis[0]) ;
396 ec->SetM20(elipAxis[1]*elipAxis[1]) ;
397 ec->SetM11(-1) ; //not yet implemented
399 TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
400 arrayTrackMatched[0]= matchedTrack[iClust];
401 ec->AddTracksMatched(arrayTrackMatched);
403 TArrayI arrayParents(parentMult,parentList);
404 ec->AddLabels(arrayParents);
407 // add the cluster to the esd object
408 esd->AddCaloCluster(ec);
412 } // cycle on clusters
414 delete [] matchedTrack;
416 esd->SetNumberOfEMCALClusters(nClustersNew);
417 //if(nClustersNew != nClusters)
418 //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
420 //Fill ESDCaloCluster with PID weights
421 AliEMCALPID *pid = new AliEMCALPID;
422 //pid->SetPrintInfo(kTRUE);
423 pid->SetReconstructor(kTRUE);
429 printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew);
433 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
435 // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
437 // Clear previous digits
438 fgDigitsArr->Delete();
441 // Read the digits from the input tree
442 TBranch *branch = digitsTree->GetBranch("EMCAL");
444 AliError("can't get the branch with the EMCAL digits !");
447 fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
448 branch->SetAddress(&fgDigitsArr);