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