AliPhysics  8b695ca (8b695ca)
AliEmcalCorrectionCellTimeCalib.cxx
Go to the documentation of this file.
1 // AliEmcalCorrectionCellTimeCalib
2 //
3 
4 #include <memory>
5 
6 #include <TObjArray.h>
7 #include <TFile.h>
8 #include "AliEMCALGeometry.h"
9 #include "AliOADBContainer.h"
10 #include "AliEMCALRecoUtils.h"
11 #include "AliAODEvent.h"
12 #include "AliDataFile.h"
13 
15 
19 
20 // Actually registers the class with the base class
22 
27  AliEmcalCorrectionComponent("AliEmcalCorrectionCellTimeCalib")
28  ,fCellTimeDistBefore(0)
29  ,fCellTimeDistAfter(0)
30  ,fCalibrateTime(kFALSE)
31  ,fCalibrateTimeL1Phase(kFALSE)
32  ,fUseAutomaticTimeCalib(1)
33 {
34 }
35 
40 {
41 }
42 
47 {
48  // Initialization
50 
51  AliWarning("Init EMCAL time calibration");
52 
53  fCalibrateTime = kTRUE;
54 
55  // init reco utils
56  if (!fRecoUtils)
58 
60 
61  return kTRUE;
62 }
63 
68 {
70 
71  if (fCreateHisto){
72  fCellTimeDistBefore = new TH1F("hCellTimeDistBefore","hCellTimeDistBefore;t_{cell} (s)",1000,-10e-6,10e-6);
74  fCellTimeDistAfter = new TH1F("hCellTimeDistAfter","hCellTimeDistAfter;t_{cell} (s)",1000,-10e-6,10e-6);
76  }
77 }
78 
83 {
85 
86  if (!fEventManager.InputEvent()) {
87  AliError("Event ptr = 0, returning");
88  return kFALSE;
89  }
90 
92 
93  // CONFIGURE THE RECO UTILS -------------------------------------------------
94 
95  // allows time calibration
96  if (fCalibrateTime)
98  else
100 
101  // allows time calibration with L1 phase
104  else
106 
107  // START PROCESSING ---------------------------------------------------------
108  // Test if cells present
109  if (fCaloCells->GetNumberOfCells()<=0)
110  {
111  AliWarning(Form("Number of EMCAL cells = %d, returning", fCaloCells->GetNumberOfCells()));
112  return kFALSE;
113  }
114 
115  // mark the cells not recalibrated
117 
118  if(fCreateHisto)
119  FillCellQA(fCellTimeDistBefore); // "before" QA
120 
121  // CELL RECALIBRATION -------------------------------------------------------
122  // cell objects will be updated
123  UpdateCells();
124 
125  if(fCreateHisto)
126  FillCellQA(fCellTimeDistAfter); // "after" QA
127 
128  return kTRUE;
129 }
130 
135 {
136  // Initialising bad channel maps
137 
138  if (!fEventManager.InputEvent())
139  return 0;
140 
141  AliInfo("Initialising time calibration map");
142 
143  // init default maps first
146 
147  Int_t runBC = fEventManager.InputEvent()->GetRunNumber();
148 
149  std::unique_ptr<AliOADBContainer> contTimeCalib;
150  std::unique_ptr<TFile> timeCalibFile;
151  if (fBasePath!="")
152  { //if fBasePath specified in the ->SetBasePath()
153  AliInfo(Form("Loading time calibration OADB from given path %s",fBasePath.Data()));
154 
155  timeCalibFile = std::unique_ptr<TFile>(TFile::Open(Form("%s/EMCALTimeCalib.root",fBasePath.Data()),"read"));
156  if (!timeCalibFile || timeCalibFile->IsZombie())
157  {
158  AliFatal(Form("EMCALTimeCalib.root was not found in the path provided: %s",fBasePath.Data()));
159  return 0;
160  }
161 
162  contTimeCalib = std::unique_ptr<AliOADBContainer>(static_cast<AliOADBContainer *>(timeCalibFile->Get("AliEMCALTimeCalib")));
163  }
164  else
165  { // Else choose the one in the $ALICE_PHYSICS directory
166  AliInfo("Loading time calibration OADB from $ALICE_PHYSICS/OADB/EMCAL");
167 
168  timeCalibFile = std::unique_ptr<TFile>(TFile::Open(AliDataFile::GetFileNameOADB("EMCAL/EMCALTimeCalib.root").data(),"read"));
169  if (!timeCalibFile || timeCalibFile->IsZombie())
170  {
171  AliFatal("OADB/EMCAL/EMCALTimeCalib.root was not found");
172  return 0;
173  }
174 
175  contTimeCalib = std::unique_ptr<AliOADBContainer>(static_cast<AliOADBContainer *>(timeCalibFile->Get("AliEMCALTimeCalib")));
176  }
177  if(!contTimeCalib){
178  AliError("No OADB container found");
179  return 0;
180  }
181  contTimeCalib->SetOwner(true);
182 
183  TObjArray *arrayBC=(TObjArray*)contTimeCalib->GetObject(runBC);
184  if (!arrayBC)
185  {
186  AliError(Form("No external time calibration set for run number: %d", runBC));
187  return 2;
188  }
189 
190  // The calibration object is accessed by specifying a pass
191  // For run 1, the actual pass is used (fFilepass, as determined in AliEmcalCorrectionComponent::GetPass())
192  // For run 2, the pass is always set to pass1 (as a convention)
193  // Other exceptions are hard-coded below
194 
195  TString pass = fFilepass;
196  if (fFilepass=="spc_calo") pass = "pass3";
197  if (fRun > 209121) pass = "pass1";
198 
199  TObjArray *arrayBCpass=(TObjArray*)arrayBC->FindObject(pass);
200  if (!arrayBCpass)
201  {
202  AliError(Form("No external time calibration set for: %d -%s", runBC,pass.Data()));
203  return 2;
204  }
205 
206  arrayBCpass->Print();
207 
208  for(Int_t i = 0; i < 4; i++)
209  {
211  if (h)
212  delete h;
213 
214  h = (TH1F*)arrayBCpass->FindObject(Form("hAllTimeAvBC%d",i));
215 
216  if (!h)
217  {
218  AliError(Form("Can not get hAllTimeAvBC%d",i));
219  continue;
220  }
221 
222  // Shift parameters for bc0 and bc1 in this pass
223  if ( pass=="spc_calo" && (i==0 || i==1) )
224  {
225  for(Int_t icell = 0; icell < h->GetNbinsX(); icell++)
226  h->SetBinContent(icell,h->GetBinContent(icell)-100);
227  }
228 
229  h->SetDirectory(0);
231  }
232 
233  return 1;
234 }
235 
240 {
241  // Initialising run-by-run L1 phase in time calibration maps
242 
243  if (!fEventManager.InputEvent())
244  return 0;
245 
246  AliInfo("Initialising run-by-run L1 phase in time calibration map");
247 
248  // init default maps first
251 
252  Int_t runBC = fEventManager.InputEvent()->GetRunNumber();
253 
254  std::unique_ptr<AliOADBContainer> contTimeCalib;
255  std::unique_ptr<TFile> timeFile;
256  if (fBasePath!="")
257  { //if fBasePath specified in the ->SetBasePath()
258  AliInfo(Form("Loading time calibration OADB from given path %s",fBasePath.Data()));
259 
260  timeFile = std::unique_ptr<TFile>(TFile::Open(Form("%s/EMCALTimeL1PhaseCalib.root",fBasePath.Data()),"read"));
261  if (!timeFile || timeFile->IsZombie())
262  {
263  AliFatal(Form("EMCALTimeL1PhaseCalib.root was not found in the path provided: %s",fBasePath.Data()));
264  return 0;
265  }
266 
267  contTimeCalib = std::unique_ptr<AliOADBContainer>(static_cast<AliOADBContainer *>(timeFile->Get("AliEMCALTimeL1PhaseCalib")));
268  }
269  else
270  { // Else choose the one in the $ALICE_PHYSICS directory
271  AliInfo("Loading L1 phase in time calibration OADB from OADB/EMCAL");
272 
273  timeFile = std::unique_ptr<TFile>(TFile::Open(AliDataFile::GetFileNameOADB("EMCAL/EMCALTimeL1PhaseCalib.root").data(),"read"));
274  if (!timeFile || timeFile->IsZombie())
275  {
276  AliFatal("OADB/EMCAL/EMCALTimeL1PhaseCalib.root was not found");
277  return 0;
278  }
279 
280  contTimeCalib = std::unique_ptr<AliOADBContainer>(static_cast<AliOADBContainer *>(timeFile->Get("AliEMCALTimeL1PhaseCalib")));
281  }
282  if(!contTimeCalib){
283  AliError("No OADB container found");
284  return 0;
285  }
286  contTimeCalib->SetOwner(true);
287 
288  TObjArray *arrayBC=(TObjArray*)contTimeCalib->GetObject(runBC);
289  if (!arrayBC)
290  {
291  AliError(Form("No external L1 phase in time calibration set for run number: %d", runBC));
292  return 2;
293  }
294 
295  // The calibration object is accessed by specifying a pass
296  // For run 2 (which is the only time L1-phase is implemented), the pass is always set to pass1 (as a convention)
297  // Other exceptions are hard-coded below
298 
299  TString pass = "pass1";
300 
301  if ( fFilepass=="muon_calo_pass1" && fRun > 209121 && fRun < 244284 )
302  pass = "pass0";//period LHC15a-m
303 
304  TObjArray *arrayBCpass=(TObjArray*)arrayBC->FindObject(pass);
305  if (!arrayBCpass)
306  {
307  AliError(Form("No external L1 phase in time calibration set for: %d -%s", runBC,pass.Data()));
308  return 2;
309  }
310 
311  arrayBCpass->Print();
312 
313 
315  if (h) delete h;
316 
317  h = (TH1C*)arrayBCpass->FindObject(Form("h%d",runBC));
318 
319  if (!h) {
320  AliFatal(Form("There is no calibration histogram h%d for this run",runBC));
321  }
322  h->SetDirectory(0);
324 
325  return 1;
326 }
327 
333 {
335 
336  if (runChanged) {
337 
338  // define what recalib parameters are needed for various switches
339  // this is based on implementation in AliEMCALRecoUtils
340  Bool_t needTimecalib = fCalibrateTime;
341  if(fRun>209121) fCalibrateTimeL1Phase = kTRUE;
342  Bool_t needTimecalibL1Phase = fCalibrateTime & fCalibrateTimeL1Phase;
343 
344  // init time calibration
345  if (needTimecalib && fUseAutomaticTimeCalib) {
346  Int_t initTC = InitTimeCalibration();
347  if (!initTC) {
348  AliError("InitTimeCalibration returned false, returning");
349  }
350  if (initTC==1) {
351  AliWarning("InitTimeCalib OK");
352  }
353  if (initTC > 1) {
354  AliWarning(Form("No external time calibration available: %d - %s", fEventManager.InputEvent()->GetRunNumber(), fFilepass.Data()));
355  }
356  }
357 
358  // init time calibration with L1 phase
359  if (needTimecalibL1Phase && fUseAutomaticTimeCalib) {
360  Int_t initTCL1Phase = InitTimeCalibrationL1Phase();
361  if (!initTCL1Phase) {
362  AliError("InitTimeCalibrationL1Phase returned false, returning");
363  }
364  if (initTCL1Phase==1) {
365  AliWarning("InitTimeCalibL1Phase OK");
366  }
367  if (initTCL1Phase > 1) {
368  AliWarning(Form("No external time calibration L1 phase available: %d - %s", fEventManager.InputEvent()->GetRunNumber(), fFilepass.Data()));
369  }
370  }
371  }
372  return runChanged;
373 }
TH1F * fCellTimeDistAfter
! cell energy distribution, after time calibration
TH1F * GetEMCALChannelTimeRecalibrationFactors(Int_t bc) const
TH1F * fCellTimeDistBefore
! cell energy distribution, before time calibration
void InitEMCALL1PhaseInTimeRecalibration()
TH1C * GetEMCALL1PhaseInTimeRecalibrationForAllSM() const
Bool_t fUseAutomaticTimeCalib
On by default the check in the OADB of the time recalibration.
AliVCaloCells * fCaloCells
! Pointer to CaloCells
Some utilities for cluster and cell treatment.
TObjArray * GetEMCALL1PhaseInTimeRecalibrationArray() const
AliEMCALRecoUtils * fRecoUtils
Pointer to RecoUtils.
void SwitchOnL1PhaseInTimeRecalibration()
int Int_t
Definition: External.C:63
Bool_t fCalibrateTimeL1Phase
flag cell time calibration with L1phase shift
Base class for correction components in the EMCal correction framework.
TString fBasePath
Base folder path to get root files.
void SetPositionAlgorithm(Int_t alg)
void SwitchOffL1PhaseInTimeRecalibration()
TList * fOutput
! List of output histograms
Bool_t fCreateHisto
Flag to make some basic histograms.
void InitEMCALTimeRecalibrationFactors()
Time calibration correction component in the EMCal correction framework.
void SetEMCALL1PhaseInTimeRecalibrationForAllSM(const TObjArray *map)
AliEmcalCorrectionEventManager fEventManager
Minimal task which inherits from AliAnalysisTaskSE and manages access to the event.
Bool_t fCalibrateTime
flag cell time calibration
TObjArray * GetEMCALTimeRecalibrationFactorsArray() const
static RegisterCorrectionComponent< AliEmcalCorrectionCellTimeCalib > reg
void SetEMCALChannelTimeRecalibrationFactors(const TObjArray *map)
bool Bool_t
Definition: External.C:53
TString fFilepass
Input data pass number.