]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Jet trigger staf which was used at EMCAL CD2 simulation;
[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
31 // --- Standard library ---
32
33 // --- AliRoot header files ---
34 #include "AliEMCALReconstructor.h"
35
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"
53 #include "AliEMCAL.h"
54 #include "AliEMCALHistoUtilities.h"
55 #include "AliESDVZERO.h"
56
57 #include "AliRunLoader.h"
58 #include "AliRun.h"
59
60 ClassImp(AliEMCALReconstructor) 
61
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
65
66 //____________________________________________________________________________
67 AliEMCALReconstructor::AliEMCALReconstructor() 
68   : fDebug(kFALSE), fList(0), fGeom(0) 
69 {
70   // ctor
71   InitRecParam();
72
73   fgRawUtils = new AliEMCALRawUtils;
74
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();
80   } else {
81     AliInfo(Form("Using default geometry in reconstruction"));
82     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
83   }
84
85   if(!fGeom) AliFatal(Form("Could not get geometry!"));
86
87
88
89 //____________________________________________________________________________
90 AliEMCALReconstructor::AliEMCALReconstructor(const AliEMCALReconstructor & rec)
91   : AliReconstructor(rec),
92     fDebug(rec.fDebug),
93     fList(rec.fList),
94     fGeom(rec.fGeom)
95 {
96   //copy ctor
97 }
98
99 //____________________________________________________________________________
100 AliEMCALReconstructor::~AliEMCALReconstructor()
101 {
102   // dtor
103   delete fGeom;
104   AliCodeTimer::Instance()->Print();
105
106
107 //____________________________________________________________________________
108 void AliEMCALReconstructor::Init()
109 {
110   // Trigger hists - Oct 24, 2007
111   fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
112 }
113
114 //____________________________________________________________________________
115 void AliEMCALReconstructor::InitRecParam() const
116 {
117   // Check if the instance of AliEMCALRecParam exists, 
118   // if not, get it from OCDB if available, otherwise create a default one
119
120  if (!fgkRecParam  && (AliCDBManager::Instance()->IsDefaultStorageSet())) {
121     AliCDBEntry *entry = (AliCDBEntry*) 
122       AliCDBManager::Instance()->Get("EMCAL/Config/RecParam");
123     if (entry) fgkRecParam =  (AliEMCALRecParam*) entry->GetObject();
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 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
134 {
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.
140
141   AliCodeTimerAuto("")
142
143   ReadDigitsArrayFromTree(digitsTree);
144
145   if(fgDigitsArr && fgDigitsArr->GetEntries()) {
146
147     AliEMCALClusterizerv1 clu;
148     clu.SetInput(digitsTree);
149     clu.SetOutput(clustersTree);
150     if(Debug())
151       clu.Digits2Clusters("deb all") ;
152     else
153       clu.Digits2Clusters("") ;
154
155   }
156 }
157
158 //____________________________________________________________________________
159 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
160
161 {
162   // Conversion from raw data to
163   // EMCAL digits.
164   // Works on a single-event basis
165
166   rawReader->Reset() ; 
167
168   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",100);
169   Int_t bufsize = 32000;
170   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
171
172   //must be done here because, in constructor, option is not yet known
173   fgRawUtils->SetOption(GetOption());
174
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());
180
181   fgRawUtils->Raw2Digits(rawReader,digitsArr);
182
183   digitsTree->Fill();
184   digitsArr->Delete();
185   delete digitsArr;
186
187 }
188
189
190 //____________________________________________________________________________
191 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
192                                     AliESDEvent* esd) const
193 {
194   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
195   // and V0 
196   // Works on the current event
197   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
198   //return;
199   const double timeScale = 1.e+11; // transition constant from sec to 0.01 ns 
200
201   //######################################################
202   //#########Calculate trigger and set trigger info###########
203   //######################################################
204  
205   AliEMCALTrigger tr;
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 
211   // at ESD
212   AliESDVZERO* vZero = esd->GetVZEROData();
213   if(vZero) {
214     tr.SetVZER0Multiplicity(vZero->GetMTotV0A() + vZero->GetMTotV0C());
215   }
216   //
217   tr.Trigger();
218
219   Float_t maxAmp2x2  = tr.Get2x2MaxAmplitude();
220   Float_t maxAmpnxn  = tr.GetnxnMaxAmplitude();
221   Float_t ampOutOfPatch2x2  = tr.Get2x2AmpOutOfPatch() ;
222   Float_t ampOutOfPatchnxn  = tr.GetnxnAmpOutOfPatch() ;
223
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();
230
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));
233
234   TVector3    pos2x2(-1,-1,-1);
235   TVector3    posnxn(-1,-1,-1);
236
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);
242   
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]);
252
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]);
260   //printf("\n");
261   tr.Print("");
262
263   esd->AddEMCALTriggerPosition(triggerPosition);
264   esd->AddEMCALTriggerAmplitudes(triggerAmplitudes);
265   // Fill trigger hists
266   AliEMCALHistoUtilities::FillTriggersListOfHists(fList,&triggerPosition,&triggerAmplitudes);
267
268   //########################################
269   //##############Fill CaloCells###############
270   //########################################
271
272   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
273   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
274   if (!branchdig) { 
275     AliError("can't get the branch with the PHOS digits !");
276     return;
277   }
278   branchdig->SetAddress(&digits);
279   digitsTree->GetEvent(0);
280   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
281   AliDebug(1,Form("%d digits",nDigits));
282
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());   
290       idignew++;
291     }
292   }
293   emcCells.SetNumberOfCells(idignew);
294   emcCells.Sort();
295
296   //------------------------------------------------------------
297   //-----------------CLUSTERS-----------------------------
298   //------------------------------------------------------------
299   TObjArray *clusters = new TObjArray(100);
300   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
301   branch->SetAddress(&clusters);
302   clustersTree->GetEvent(0);
303
304   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
305   AliDebug(1,Form("%d clusters",nClusters));
306   esd->SetFirstEMCALCluster(esd->GetNumberOfCaloClusters()); // Put after Phos clusters 
307
308
309   //######################################################
310   //#######################TRACK MATCHING###############
311   //######################################################
312   //Fill list of integers, each one is index of track to which the cluster belongs.
313
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
318   
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;
326   } 
327   
328   //########################################
329   //##############Fill CaloClusters#############
330   //########################################
331
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
337     Float_t xyz[3];
338     TVector3 gpos;
339     clust->GetGlobalPosition(gpos);
340     for (Int_t ixyz=0; ixyz<3; ixyz++) 
341       xyz[ixyz] = gpos[ixyz];
342     
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();
350     Float_t elipAxis[2];
351     clust->GetElipsAxis(elipAxis);
352     
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);
362         else
363           timeList[newdigitMult] = 65535;
364         digiList[newdigitMult] = (UShort_t)(digitInts[iDigit]);
365         newdigitMult++;
366       }
367       else if (clust->GetClusterType() != AliESDCaloCluster::kEMCALPseudoCluster)
368         Warning("FillESD()","Negative or 0 digit amplitude in cluster");
369     }
370     
371     if(newdigitMult > 0) { // accept cluster if it has some digit
372       nClustersNew++;
373       
374       amplList.Set(newdigitMult);
375       timeList.Set(newdigitMult);
376       digiList.Set(newdigitMult);
377
378       //Primaries
379       Int_t  parentMult  = 0;
380       Int_t *parentList =  clust->GetParents(parentMult);
381     
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);
390     
391       if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1){
392
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
398
399        TArrayI arrayTrackMatched(1);// Only one track, temporal solution. 
400        arrayTrackMatched[0]= matchedTrack[iClust]; 
401        ec->AddTracksMatched(arrayTrackMatched); 
402          
403        TArrayI arrayParents(parentMult,parentList); 
404        ec->AddLabels(arrayParents);
405       } 
406       
407       // add the cluster to the esd object
408       esd->AddCaloCluster(ec);
409       delete ec;
410     }
411
412   } // cycle on clusters
413
414   delete [] matchedTrack;
415
416   esd->SetNumberOfEMCALClusters(nClustersNew);
417   //if(nClustersNew != nClusters) 
418   //printf(" ##### nClusters %i -> new %i ##### \n", nClusters, nClustersNew );
419   
420   //Fill ESDCaloCluster with PID weights
421   AliEMCALPID *pid = new AliEMCALPID;
422   //pid->SetPrintInfo(kTRUE);
423   pid->SetReconstructor(kTRUE);
424   pid->RunPID(esd);
425
426   delete pid;
427   delete clusters;
428
429   printf(" ## AliEMCALReconstructor::FillESD() is ended : ncl %i -> %i ### \n ",nClusters, nClustersNew); 
430   //  assert(0);
431 }
432
433 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
434 {
435   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
436   if(fgDigitsArr) {
437     // Clear previous digits 
438     fgDigitsArr->Delete();
439     delete fgDigitsArr;
440   }
441   // Read the digits from the input tree
442   TBranch *branch = digitsTree->GetBranch("EMCAL");
443   if (!branch) { 
444     AliError("can't get the branch with the EMCAL digits !");
445     return;
446   }
447   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
448   branch->SetAddress(&fgDigitsArr);
449   branch->GetEntry(0);
450 }
451