]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
increase array size to avoid seg faults in testbeam reconstruction
[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 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  && (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
134 //____________________________________________________________________________
135 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
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.
141   // Works on the current event.
142
143   AliCodeTimerAuto("")
144
145   ReadDigitsArrayFromTree(digitsTree);
146
147   fgClusterizer->SetOutput(clustersTree);
148
149   if(fgDigitsArr && fgDigitsArr->GetEntries()) {
150
151     fgClusterizer->SetInput(digitsTree);
152     
153     if(Debug())
154       fgClusterizer->Digits2Clusters("deb all") ;
155     else
156       fgClusterizer->Digits2Clusters("");
157     
158     fgClusterizer->Clear();
159
160   }
161
162 }
163
164 //____________________________________________________________________________
165 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
166
167 {
168   // Conversion from raw data to
169   // EMCAL digits.
170   // Works on a single-event basis
171
172   rawReader->Reset() ; 
173
174   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
175   Int_t bufsize = 32000;
176   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
177
178   //must be done here because, in constructor, option is not yet known
179   fgRawUtils->SetOption(GetOption());
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
187   fgRawUtils->Raw2Digits(rawReader,digitsArr);
188
189   digitsTree->Fill();
190   digitsArr->Delete();
191   delete digitsArr;
192
193 }
194
195
196 //____________________________________________________________________________
197 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
198                                     AliESDEvent* esd) const
199 {
200   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
201   // and V0 
202   // Works on the current event
203   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
204   //return;
205
206   //######################################################
207   //#########Calculate trigger and set trigger info###########
208   //######################################################
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   //
222   tr.Trigger();
223
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
229   Int_t iSM2x2      = tr.Get2x2SuperModule();
230   Int_t iSMnxn      = tr.GetnxnSuperModule();
231   Int_t iModulePhi2x2 = tr.Get2x2ModulePhi();
232   Int_t iModulePhinxn = tr.GetnxnModulePhi();
233   Int_t iModuleEta2x2 = tr.Get2x2ModuleEta();
234   Int_t iModuleEtanxn = tr.GetnxnModuleEta();
235
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));
238
239   TVector3    pos2x2(-1,-1,-1);
240   TVector3    posnxn(-1,-1,-1);
241
242   Int_t iAbsId2x2 = fGeom->GetAbsCellIdFromCellIndexes( iSM2x2, iModulePhi2x2, iModuleEta2x2) ; // should be changed to Module
243   Int_t iAbsIdnxn = fGeom->GetAbsCellIdFromCellIndexes( iSMnxn, iModulePhinxn, iModuleEtanxn) ;
244   fGeom->GetGlobal(iAbsId2x2, pos2x2);
245   fGeom->GetGlobal(iAbsIdnxn, posnxn);
246   //printf(" iAbsId2x2 %i iAbsIdnxn %i \n", iAbsId2x2, iAbsIdnxn);
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) ;   
254   triggerPosition[5] = posnxn(2) ;
255   //printf(" triggerPosition ");
256   //for(int i=0; i<6; i++) printf(" %i %f : ", i, triggerPosition[i]);
257
258   TArrayF triggerAmplitudes(4);
259   triggerAmplitudes[0] = maxAmp2x2 ;   
260   triggerAmplitudes[1] = ampOutOfPatch2x2 ;    
261   triggerAmplitudes[2] = maxAmpnxn ;   
262   triggerAmplitudes[3] = ampOutOfPatchnxn ;   
263   //printf("\n triggerAmplitudes ");
264   //for(int i=0; i<4; i++) printf(" %i %f : ", i, triggerAmplitudes[i]);
265   //printf("\n");
266   //tr.Print("");
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   }
284   esd->AddEMCALTriggerPosition(triggerPosition);
285   esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
286   // Fill trigger hists
287   AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
288
289   //########################################
290   //##############Fill CaloCells###############
291   //########################################
292
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
329
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();
346     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
347   } 
348   
349   //########################################
350   //##############Fill CaloClusters#############
351   //########################################
352   esd->SetNumberOfEMCALClusters(nClusters);
353   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
354     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
355     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
356     if (Debug()) clust->Print();
357     // Get information from EMCAL reconstruction points
358     Float_t xyz[3];
359     TVector3 gpos;
360     clust->GetGlobalPosition(gpos);
361     for (Int_t ixyz=0; ixyz<3; ixyz++)
362       xyz[ixyz] = gpos[ixyz];
363     Float_t elipAxis[2];
364     clust->GetElipsAxis(elipAxis);
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++;
381       }
382     }
383     absIdList.Set(newCellMult);
384     //Uncomment when unfolding is done
385     //fracList.Set(newCellMult);
386
387     if(newCellMult > 0) { // accept cluster if it has some digit
388       nClustersNew++;
389       //Primaries
390       Int_t  parentMult  = 0;
391       Int_t *parentList =  clust->GetParents(parentMult);
392       // fills the ESDCaloCluster
393       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
394       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
395       ec->SetPosition(xyz);
396       ec->SetE(clust->GetEnergy());
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]) ;
414       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
415       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
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
423       // add the cluster to the esd object
424      esd->AddCaloCluster(ec);
425       delete ec;
426       delete [] newAbsIdList ;
427       //delete [] newFracList ;
428    }
429  } // cycle on clusters
430
431  delete [] matchedTrack;
432
433  esd->SetNumberOfEMCALClusters(nClustersNew);
434  //if(nClustersNew != nClusters)
435  //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
436
437  //Fill ESDCaloCluster with PID weights
438   AliEMCALPID *pid = new AliEMCALPID;
439   //pid->SetPrintInfo(kTRUE);
440   pid->SetReconstructor(kTRUE);
441   pid->RunPID(esd);
442   delete pid;
443
444   delete digits;
445   delete clusters;
446
447   // printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew); 
448 }
449
450 //__________________________________________________________________________
451 void 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 }
469