]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/AliEMCALReconstructor.cxx
Skip Calibration events by default. If not skip do not clusterize, just fill the...
[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   }
245   else{
246         AliDebug(1," Calibration Event, skip!");
247         printf("**** AliEMCALReconstructor::ConvertDigits() Calibration Event, skip!!\n");
248   }
249         
250   digitsTree->Fill();
251   digitsArr->Delete();
252   digitsTrg->Delete();
253   delete digitsArr;
254   delete digitsTrg;
255
256 }
257
258
259 //____________________________________________________________________________
260 void AliEMCALReconstructor::FillESD(TTree* digitsTree, TTree* clustersTree, 
261                                     AliESDEvent* esd) const
262 {
263   // Called by AliReconstruct after Reconstruct() and global tracking and vertexing 
264   // and V0 
265   // Works on the current event
266   //  printf(" ## AliEMCALReconstructor::FillESD() is started ### \n ");
267   //return;
268
269   //########################################
270   //##############Fill CaloCells###############
271   //########################################
272
273   TClonesArray *digits = new TClonesArray("AliEMCALDigit",1000);
274   TBranch *branchdig = digitsTree->GetBranch("EMCAL");
275   if (!branchdig) { 
276     AliError("can't get the branch with the EMCAL digits !");
277     return;
278   }
279   branchdig->SetAddress(&digits);
280   digitsTree->GetEvent(0);
281   Int_t nDigits = digits->GetEntries(), idignew = 0 ;
282   AliDebug(1,Form("%d digits",nDigits));
283
284   AliESDCaloCells &emcCells = *(esd->GetEMCALCells());
285   emcCells.CreateContainer(nDigits);
286   emcCells.SetType(AliESDCaloCells::kEMCALCell);
287   Float_t energy = 0;
288   for (Int_t idig = 0 ; idig < nDigits ; idig++) {
289     const AliEMCALDigit * dig = (const AliEMCALDigit*)digits->At(idig);
290     if(dig->GetAmplitude() > 0 ){
291           energy = (static_cast<AliEMCALClusterizerv1*> (fgClusterizer))->Calibrate(dig->GetAmplitude(),dig->GetTime(),dig->GetId()); //TimeR or Time?
292           if(energy > 0){ //Digits tagged as bad (dead, hot, not alive) are set to 0 in calibrate, remove them  
293                   emcCells.SetCell(idignew,dig->GetId(),energy, dig->GetTime());   
294                   idignew++;
295           }
296     }
297   }
298   emcCells.SetNumberOfCells(idignew);
299   emcCells.Sort();
300
301   //------------------------------------------------------------
302   //-----------------CLUSTERS-----------------------------
303   //------------------------------------------------------------
304   clustersTree->SetBranchStatus("*",0); //disable all branches
305   clustersTree->SetBranchStatus("EMCALECARP",1); //Enable only the branch we need
306
307   TObjArray *clusters = new TObjArray(100);
308   TBranch *branch = clustersTree->GetBranch("EMCALECARP");
309   branch->SetAddress(&clusters);
310   branch->GetEntry(0);
311   //clustersTree->GetEvent(0);
312
313   Int_t nClusters = clusters->GetEntries(),  nClustersNew=0;
314   AliDebug(1,Form("%d clusters",nClusters));
315
316   //######################################################
317   //#######################TRACK MATCHING###############
318   //######################################################
319   //Fill list of integers, each one is index of track to which the cluster belongs.
320
321   // step 1 - initialize array of matched track indexes
322   Int_t *matchedTrack = new Int_t[nClusters];
323   for (Int_t iclus = 0; iclus < nClusters; iclus++)
324     matchedTrack[iclus] = -1;  // neg. index --> no matched track
325   
326   // step 2, change the flag for all matched clusters found in tracks
327   Int_t iemcalMatch = -1;
328   Int_t endtpc = esd->GetNumberOfTracks();
329   for (Int_t itrack = 0; itrack < endtpc; itrack++) {
330     AliESDtrack * track = esd->GetTrack(itrack) ; // retrieve track
331     iemcalMatch = track->GetEMCALcluster();
332     if(iemcalMatch >= 0) matchedTrack[iemcalMatch] = itrack;
333   } 
334   
335   //########################################
336   //##############Fill CaloClusters#############
337   //########################################
338   for (Int_t iClust = 0 ; iClust < nClusters ; iClust++) {
339     const AliEMCALRecPoint * clust = (const AliEMCALRecPoint*)clusters->At(iClust);
340     //if(clust->GetClusterType()== AliESDCaloCluster::kEMCALClusterv1) nRP++; else nPC++;
341     if (Debug()) clust->Print();
342     // Get information from EMCAL reconstruction points
343     Float_t xyz[3];
344     TVector3 gpos;
345     clust->GetGlobalPosition(gpos);
346     for (Int_t ixyz=0; ixyz<3; ixyz++)
347       xyz[ixyz] = gpos[ixyz];
348     Float_t elipAxis[2];
349     clust->GetElipsAxis(elipAxis);
350        //Create digits lists
351     Int_t cellMult = clust->GetMultiplicity();
352     //TArrayS digiList(digitMult);
353     Float_t *amplFloat = clust->GetEnergiesList();
354     Int_t   *digitInts = clust->GetAbsId();
355     TArrayS absIdList(cellMult);
356     TArrayD fracList(cellMult);
357
358     Int_t newCellMult = 0;
359     for (Int_t iCell=0; iCell<cellMult; iCell++) {
360       if (amplFloat[iCell] > 0) {
361       absIdList[newCellMult] = (UShort_t)(digitInts[iCell]);
362       //Uncomment when unfolding is done
363       //if(emcCells.GetCellAmplitude(digitInts[iCell])>0)
364       //fracList[newCellMult] = amplFloat[iCell]/(emcCells.GetCellAmplitude(digitInts[iCell])*calibration);//get cell calibration value 
365       //else
366       fracList[newCellMult] = 0; 
367       newCellMult++;
368       }
369     }
370
371     absIdList.Set(newCellMult);
372     fracList.Set(newCellMult);
373     
374     if(newCellMult > 0) { // accept cluster if it has some digit
375       nClustersNew++;
376       //Primaries
377       Int_t  parentMult  = 0;
378       Int_t *parentList =  clust->GetParents(parentMult);
379       // fills the ESDCaloCluster
380       AliESDCaloCluster * ec = new AliESDCaloCluster() ;
381       ec->SetClusterType(AliESDCaloCluster::kEMCALClusterv1);
382       ec->SetPosition(xyz);
383       ec->SetE(clust->GetEnergy());
384                 
385           //Distance to the nearest bad crystal
386           ec->SetDistanceToBadChannel(clust->GetDistanceToBadTower()); 
387
388       ec->SetNCells(newCellMult);
389       //Change type of list from short to ushort
390       UShort_t *newAbsIdList  = new UShort_t[newCellMult];
391       Double_t *newFracList  = new Double_t[newCellMult];
392       for(Int_t i = 0; i < newCellMult ; i++) {
393         newAbsIdList[i]=absIdList[i];
394         newFracList[i]=fracList[i];
395       }
396       ec->SetCellsAbsId(newAbsIdList);
397       ec->SetCellsAmplitudeFraction(newFracList);
398       ec->SetClusterDisp(clust->GetDispersion());
399       ec->SetClusterChi2(-1); //not yet implemented
400       ec->SetM02(elipAxis[0]*elipAxis[0]) ;
401       ec->SetM20(elipAxis[1]*elipAxis[1]) ;
402       ec->SetTOF(clust->GetTime()) ; //time-of-fligh
403       ec->SetNExMax(clust->GetNExMax());          //number of local maxima
404       TArrayI arrayTrackMatched(1);// Only one track, temporal solution.
405       arrayTrackMatched[0]= matchedTrack[iClust];
406       ec->AddTracksMatched(arrayTrackMatched);
407
408       TArrayI arrayParents(parentMult,parentList);
409       ec->AddLabels(arrayParents);
410
411       // add the cluster to the esd object
412       esd->AddCaloCluster(ec);
413       delete ec;
414       delete [] newAbsIdList ;
415       delete [] newFracList ;
416    }
417  } // cycle on clusters
418
419  delete [] matchedTrack;
420
421  //Fill ESDCaloCluster with PID weights
422  AliEMCALPID *pid = new AliEMCALPID;
423  //pid->SetPrintInfo(kTRUE);
424  pid->SetReconstructor(kTRUE);
425  pid->RunPID(esd);
426  delete pid;
427   
428  delete digits;
429  delete clusters;
430   
431  //Store EMCAL misalignment matrixes
432  FillMisalMatrixes(esd) ;
433
434 }
435
436 //==================================================================================
437 void AliEMCALReconstructor::FillMisalMatrixes(AliESDEvent* esd)const{
438         //Store EMCAL matrixes in ESD Header
439         
440         //Check, if matrixes was already stored
441         for(Int_t sm = 0 ; sm < fGeom->GetNumberOfSuperModules(); sm++){
442                 if(esd->GetEMCALMatrix(sm)!=0)
443                         return ;
444         }
445         
446         //Create and store matrixes
447         if(!gGeoManager){
448                 AliError("Can not store misal. matrixes: no gGeoManager! \n") ;
449                 return ;
450         }
451         //Note, that owner of copied marixes will be header
452         char path[255] ;
453         TGeoHMatrix * m = 0x0;
454         for(Int_t sm = 0; sm < fGeom->GetNumberOfSuperModules(); sm++){
455                 sprintf(path,"/ALIC_1/XEN1_1/SMOD_%d",sm+1) ; //In Geometry modules numbered 1,2,.,5
456                 if(sm >= 10) sprintf(path,"/ALIC_1/XEN1_1/SM10_%d",sm-10+1) ;
457                 
458                 if (gGeoManager->CheckPath(path)){
459                         gGeoManager->cd(path);
460                         m = gGeoManager->GetCurrentMatrix() ;
461 //                      printf("================================================= \n");
462 //                      printf("AliEMCALReconstructor::FixMisalMatrixes(), sm %d, \n",sm);
463 //                      m->Print("");
464                         esd->SetEMCALMatrix(new TGeoHMatrix(*m),sm) ;
465 //                      printf("================================================= \n");
466                 }
467                 else{
468                         esd->SetEMCALMatrix(NULL,sm) ;
469                 }
470         }
471 }
472
473
474
475 //__________________________________________________________________________
476 void AliEMCALReconstructor::ReadDigitsArrayFromTree(TTree *digitsTree) const
477 {
478   // See AliEMCALClusterizer::SetInput(TTree *digitsTree);
479   if(fgDigitsArr) {
480     // Clear previous digits 
481     fgDigitsArr->Delete();
482     delete fgDigitsArr;
483   }
484   // Read the digits from the input tree
485   TBranch *branch = digitsTree->GetBranch("EMCAL");
486   if (!branch) { 
487     AliError("can't get the branch with the EMCAL digits !");
488     return;
489   }
490   fgDigitsArr = new TClonesArray("AliEMCALDigit",100);
491   branch->SetAddress(&fgDigitsArr);
492   branch->GetEntry(0);
493 }
494
495