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