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