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