]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Changes for report #69974: Virtual class for calorimeter analysis objects
[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 "AliEMCALClusterizerNxN.h"
49 #include "AliEMCALRecPoint.h"
50 #include "AliEMCALPID.h"
51 #include "AliEMCALTrigger.h"
52 #include "AliRawReader.h"
53 #include "AliCDBEntry.h"
54 #include "AliCDBManager.h"
55 #include "AliEMCALGeometry.h"
56 #include "AliEMCAL.h"
57 #include "AliESDVZERO.h"
58 #include "AliCDBManager.h"
59 #include "AliRunLoader.h"
60 #include "AliRun.h"
61 #include "AliEMCALTriggerData.h"
62 #include "AliEMCALTriggerElectronics.h"
63 #include "AliEMCALTriggerDCSConfigDB.h"
64 #include "AliEMCALTriggerDCSConfig.h"
65
66 ClassImp(AliEMCALReconstructor) 
67
68 const AliEMCALRecParam* AliEMCALReconstructor::fgkRecParam = 0;  // EMCAL rec. parameters
69 AliEMCALRawUtils* AliEMCALReconstructor::fgRawUtils = 0;   // EMCAL raw utilities class
70 AliEMCALClusterizer* AliEMCALReconstructor::fgClusterizer = 0;   // EMCAL clusterizer class
71 TClonesArray*     AliEMCALReconstructor::fgDigitsArr = 0;  // shoud read just once at event
72 AliEMCALTriggerElectronics* AliEMCALReconstructor::fgTriggerProcessor = 0x0;
73 //____________________________________________________________________________
74 AliEMCALReconstructor::AliEMCALReconstructor() 
75   : fDebug(kFALSE), fList(0), fGeom(0),fCalibData(0),fPedestalData(0) 
76 {
77   // ctor
78
79   fgRawUtils = new AliEMCALRawUtils;
80
81   //To make sure we match with the geometry in a simulation file,
82   //let's try to get it first.  If not, take the default geometry
83   AliRunLoader *rl = AliRunLoader::Instance();
84   if (rl->GetAliRun() && rl->GetAliRun()->GetDetector("EMCAL")) {
85     fGeom = dynamic_cast<AliEMCAL*>(rl->GetAliRun()->GetDetector("EMCAL"))->GetGeometry();
86   } else {
87     AliInfo(Form("Using default geometry in reconstruction"));
88     fGeom =  AliEMCALGeometry::GetInstance(AliEMCALGeometry::GetDefaultGeometryName());
89   }
90
91   //Get calibration parameters  
92   if(!fCalibData)
93     {
94                 AliCDBEntry *entry = (AliCDBEntry*) 
95                 AliCDBManager::Instance()->Get("EMCAL/Calib/Data");
96                 if (entry) fCalibData =  (AliEMCALCalibData*) entry->GetObject();
97     }
98         
99   if(!fCalibData)
100                 AliFatal("Calibration parameters not found in CDB!");
101         
102   //Get calibration parameters  
103   if(!fPedestalData)
104     {
105                 AliCDBEntry *entry = (AliCDBEntry*) 
106                 AliCDBManager::Instance()->Get("EMCAL/Calib/Pedestals");
107                 if (entry) fPedestalData =  (AliCaloCalibPedestal*) entry->GetObject();
108     }
109         
110   if(!fPedestalData)
111     AliFatal("Dead map not found in CDB!");
112         
113   InitClusterizer();
114                 
115   if(!fGeom) AliFatal(Form("Could not get geometry!"));
116
117   AliEMCALTriggerDCSConfigDB* dcsConfigDB = AliEMCALTriggerDCSConfigDB::Instance();
118
119   const AliEMCALTriggerDCSConfig* dcsConfig = dcsConfigDB->GetTriggerDCSConfig();
120
121   if (!dcsConfig) AliFatal("No Trigger DCS Configuration from OCDB!");
122   fgTriggerProcessor = new AliEMCALTriggerElectronics( dcsConfig );
123
124
125 //____________________________________________________________________________
126 AliEMCALReconstructor::~AliEMCALReconstructor()
127 {
128   // dtor
129   delete fGeom;
130   delete fgRawUtils;
131   delete fgClusterizer;
132   delete fgTriggerProcessor;
133
134   AliCodeTimer::Instance()->Print();
135
136
137 // //____________________________________________________________________________
138 // void AliEMCALReconstructor::Init()
139 // {
140 //   // Trigger hists - Oct 24, 2007
141 //    fList = AliEMCALHistoUtilities::GetTriggersListOfHists(kTRUE);
142 // }
143
144 //____________________________________________________________________________
145 void AliEMCALReconstructor::InitClusterizer() 
146 {
147   //Init the clusterizer with geometry and calibration pointers, avoid doing it twice.
148   
149   AliEMCALRecParam *recParam = NULL;
150   AliCDBEntry *entry = (AliCDBEntry*) 
151   AliCDBManager::Instance()->Get("EMCAL/Calib/RecoParam");
152   //Get The reco param for the default event specie
153   if (entry) 
154    recParam = (AliEMCALRecParam*)((TObjArray *) entry->GetObject())->At(0);
155     
156   if(!recParam)  
157     AliFatal("RecoParam not found in CDB!");
158     
159   if (recParam->GetClusterizerFlag() == AliEMCALRecParam::kClusterizerv1)
160   {
161     fgClusterizer = new AliEMCALClusterizerv1(fGeom, fCalibData,fPedestalData); 
162   }
163   else
164   {
165     fgClusterizer = new AliEMCALClusterizerNxN(fGeom, fCalibData,fPedestalData); 
166   }
167   
168 }
169
170 //____________________________________________________________________________
171 void AliEMCALReconstructor::Reconstruct(TTree* digitsTree, TTree* clustersTree) const
172 {
173   // method called by AliReconstruction; 
174   // Only the clusterization is performed,; the rest of the reconstruction is done in FillESD because the track
175   // segment maker needs access to the AliESD object to retrieve the tracks reconstructed by 
176   // the global tracking.
177   // Works on the current event.
178
179   AliCodeTimerAuto("",0)
180
181   ReadDigitsArrayFromTree(digitsTree);
182   
183   fgClusterizer->InitParameters();
184   fgClusterizer->SetOutput(clustersTree);
185  
186   AliEMCALTriggerData* trgData = new AliEMCALTriggerData();
187         
188   Int_t bufferSize = 32000;
189
190   if (TBranch* triggerBranch = clustersTree->GetBranch("EMTRG"))
191                 triggerBranch->SetAddress(&trgData);
192         else
193                 clustersTree->Branch("EMTRG","AliEMCALTriggerData",&trgData,bufferSize);
194
195   TClonesArray *trgDigits = new TClonesArray("AliEMCALRawDigit",1000);
196   TBranch *branchdig = digitsTree->GetBranch("EMTRG");
197   if (!branchdig) 
198   { 
199           AliError("Can't get the branch with the EMCAL trigger digits !");
200           return;
201   }
202
203   branchdig->SetAddress(&trgDigits);
204   branchdig->GetEntry(0);
205         
206   //Skip clusterization of LED events
207   if (GetRecParam()->GetEventSpecie()!=AliRecoParam::kCalib){
208
209           Int_t v0M[2] = {0,0};
210           fgTriggerProcessor->Digits2Trigger(trgDigits, v0M, trgData);
211         
212  
213           if(fgDigitsArr && fgDigitsArr->GetEntries()) {
214
215                   fgClusterizer->SetInput(digitsTree);
216     
217                   if(Debug())
218                           fgClusterizer->Digits2Clusters("deb all") ;
219                   else
220                           fgClusterizer->Digits2Clusters("");
221     
222                   fgClusterizer->Clear();
223
224           }//digits array exists and has somethind
225   }//not a LED event
226         
227   clustersTree->Fill(); 
228   trgDigits->Delete();
229   delete trgDigits; trgDigits = 0x0;
230   delete trgData;   trgData   = 0x0;
231
232 }
233
234 //____________________________________________________________________________
235 void AliEMCALReconstructor::ConvertDigits(AliRawReader* rawReader, TTree* digitsTree) const
236
237 {
238   // Conversion from raw data to
239   // EMCAL digits.
240   // Works on a single-event basis
241
242   rawReader->Reset() ; 
243
244   TClonesArray *digitsArr = new TClonesArray("AliEMCALDigit",200);
245   TClonesArray *digitsTrg = new TClonesArray("AliEMCALRawDigit", 200);
246
247   Int_t bufsize = 32000;
248   digitsTree->Branch("EMCAL", &digitsArr, bufsize);
249   digitsTree->Branch("EMTRG", &digitsTrg, bufsize);
250         
251   //Skip calibration events do the rest
252   Bool_t doFit = kTRUE;
253   if ( !(GetRecParam()->FitLEDEvents()) && GetRecParam()->GetEventSpecie()==AliRecoParam::kCalib) doFit = kFALSE;
254   if (doFit){
255           //must be done here because, in constructor, option is not yet known
256           fgRawUtils->SetOption(GetOption());
257
258           fgRawUtils->SetRawFormatHighLowGainFactor(GetRecParam()->GetHighLowGainFactor());
259           fgRawUtils->SetRawFormatOrder(GetRecParam()->GetOrderParameter());
260           fgRawUtils->SetRawFormatTau(GetRecParam()->GetTau());
261           fgRawUtils->SetNoiseThreshold(GetRecParam()->GetNoiseThreshold());
262           fgRawUtils->SetNPedSamples(GetRecParam()->GetNPedSamples());
263           fgRawUtils->SetRemoveBadChannels(GetRecParam()->GetRemoveBadChannels());
264           fgRawUtils->SetFittingAlgorithm(GetRecParam()->GetFittingAlgorithm());
265           fgRawUtils->SetFALTROUsage(GetRecParam()->UseFALTRO());
266           fgRawUtils->SetTimeMin(GetRecParam()->GetTimeMin());
267           fgRawUtils->SetTimeMax(GetRecParam()->GetTimeMax());
268         
269           fgRawUtils->Raw2Digits(rawReader,digitsArr,fPedestalData,digitsTrg);
270   }//skip calibration event
271   else{
272         AliDebug(1," Calibration Event, skip!");
273   }
274         
275   digitsTree->Fill();
276   digitsArr->Delete();
277   digitsTrg->Delete();
278   delete digitsArr;
279   delete digitsTrg;
280
281 }
282
283
284 //____________________________________________________________________________
285 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
286                                     AliESDEvent* esd) const
287 {
288   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
289   // and V0 
290   // Works on the current event
291   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
292   //return;
293
294   //########################################
295   //##############Fill CaloCells###############
296   //########################################
297
298   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
299   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
300   if (!branchdig) { 
301     AliError("can't get the branch with the EMCAL digits !");
302     return;
303   }
304   branchdig->SetAddress(&digits);
305   digitsTree->GetEvent(0);
306   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
307   AliDebug(1,Form("%d digits",nDigits));
308
309   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
310   emcCells.CreateContainer(nDigits);
311   emcCells.SetType(AliVCaloCells::kEMCALCell);
312   Float_t energy = 0;
313   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
314     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
315     if(dig->GetAmplitude() > 0 ){
316           energy = fgClusterizer->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
317           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
318                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
319                   idignew++;
320           }
321     }
322   }
323   emcCells.SetNumberOfCells(idignew);
324   emcCells.Sort();
325
326   //------------------------------------------------------------
327   //-----------------CLUSTERS-----------------------------
328   //------------------------------------------------------------
329   clustersTree->SetBranchStatus("*",0); //disable all branches
330   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
331
332   TObjArray *clusters = new TObjArray(100);
333   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
334   branch->SetAddress(&clusters);
335   branch->GetEntry(0);
336   //clustersTree->GetEvent(0);
337
338   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
339   AliDebug(1,Form("%d clusters",nClusters));
340
341   //######################################################
342   //#######################TRACK MATCHING###############
343   //######################################################
344   //Fill list of integers, each one is index of track to which the cluster belongs.
345
346   // step 1 - initialize array of matched track indexes
347   Int_t *matchedTrack = new Int_t[nClusters];
348   for (Int_t iclus = 0; iclus < nClusters; iclus++)
349     matchedTrack[iclus] = -1;  // neg. index --> no matched track
350   
351   // step 2, change the flag for all matched clusters found in tracks
352   Int_t iemcalMatch = -1;
353   Int_t endtpc = esd->GetNumberOfTracks();
354   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
355     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
356     iemcalMatch = track->GetEMCALcluster();
357     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
358   } 
359   
360   //########################################
361   //##############Fill CaloClusters#############
362   //########################################
363   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
364     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
365     //if(clust->GetClusterType()== AliVCluster::kEMCALClusterv1) nRP++; else nPC++;
366     if (Debug()) clust->Print();
367     // Get information from EMCAL reconstruction points
368     Float_t xyz[3];
369     TVector3 gpos;
370     clust->GetGlobalPosition(gpos);
371     for (Int_t ixyz=0; ixyz<3; ixyz++)
372       xyz[ixyz] = gpos[ixyz];
373     Float_t elipAxis[2];
374     clust->GetElipsAxis(elipAxis);
375        //Create digits lists
376     Int_t cellMult = clust->GetMultiplicity();
377     //TArrayS digiList(digitMult);
378     Float_t *amplFloat = clust->GetEnergiesList();
379     Int_t   *digitInts = clust->GetAbsId();
380     TArrayS absIdList(cellMult);
381     TArrayD fracList(cellMult);
382
383     Int_t newCellMult = 0;
384     for (Int_t iCell=0; iCell<cellMult; iCell++) {
385       if (amplFloat[iCell] > 0) {
386       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
387       //Uncomment when unfolding is done
388       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
389       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
390       //else
391       fracList[newCellMult] = 0; 
392       newCellMult++;
393       }
394     }
395
396     absIdList.Set(newCellMult);
397     fracList.Set(newCellMult);
398     
399     if(newCellMult > 0) { // accept cluster if it has some digit
400       nClustersNew++;
401       //Primaries
402       Int_t  parentMult  = 0;
403       Int_t *parentList =  clust->GetParents(parentMult);
404       // fills the ESDCaloCluster
405       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
406       ec->SetType(AliVCluster::kEMCALClusterv1);
407       ec->SetPosition(xyz);
408       ec->SetE(clust->GetEnergy());
409                 
410       //Distance to the nearest bad crystal
411       ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
412
413       ec->SetNCells(newCellMult);
414       //Change type of list from short to ushort
415       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
416       Double_t *newFracList   = new Double_t[newCellMult];
417       for(Int_t i = 0; i < newCellMult ; i++) {
418         newAbsIdList[i]=absIdList[i];
419         newFracList[i]=fracList[i];
420       }
421       ec->SetCellsAbsId(newAbsIdList);
422       ec->SetCellsAmplitudeFraction(newFracList);
423       ec->SetDispersion(clust->GetDispersion());
424       ec->SetChi2(-1); //not yet implemented
425       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
426       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
427       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
428       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
429       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
430       arrayTrackMatched[0]= matchedTrack[iClust];
431       ec->AddTracksMatched(arrayTrackMatched);
432
433       TArrayI arrayParents(parentMult,parentList);
434       ec->AddLabels(arrayParents);
435
436       // add the cluster to the esd object
437       esd->AddCaloCluster(ec);
438       delete ec;
439       delete [] newAbsIdList ;
440       delete [] newFracList ;
441    }
442  } // cycle on clusters
443
444  delete [] matchedTrack;
445
446  //Fill ESDCaloCluster with PID weights
447  AliEMCALPID *pid = new AliEMCALPID;
448  //pid->SetPrintInfo(kTRUE);
449  pid->SetReconstructor(kTRUE);
450  pid->RunPID(esd);
451  delete pid;
452   
453  delete digits;
454  delete clusters;
455   
456  //Store EMCAL misalignment matrixes
457  FillMisalMatrixes(esd) ;
458
459 }
460
461 //==================================================================================
462 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
463         //Store EMCAL matrixes in ESD Header
464         
465         //Check, if matrixes was already stored
466         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
467                 if(esd->GetEMCALMatrix(sm)!=0)
468                         return ;
469         }
470         
471         //Create and store matrixes
472         if(!gGeoManager){
473                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
474                 return ;
475         }
476         //Note, that owner of copied marixes will be header
477         char path[255] ;
478         TGeoHMatrix * m = 0x0;
479         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
480                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
481                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
482                 
483                 if (gGeoManager->CheckPath(path)){
484                         gGeoManager->cd(path);
485                         m = gGeoManager->GetCurrentMatrix() ;
486 //                      printf("================================================= \n");
487 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
488 //                      m->Print("");
489                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
490 //                      printf("================================================= \n");
491                 }
492                 else{
493                         esd->SetEMCALMatrix(NULL,sm) ;
494                 }
495         }
496 }
497
498
499
500 //__________________________________________________________________________
501 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
502 {
503   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
504   if(fgDigitsArr) {
505     // Clear previous digits 
506     fgDigitsArr->Delete();
507     delete fgDigitsArr;
508   }
509   // Read the digits from the input tree
510   TBranch *branch = digitsTree->GetBranch("EMCAL");
511   if (!branch) { 
512     AliError("can't get the branch with the EMCAL digits !");
513     return;
514   }
515   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
516   branch->SetAddress(&fgDigitsArr);
517   branch->GetEntry(0);
518 }
519
520