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