Belle II Software  release-08-01-10
TOPLaserCalibratorModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 // Own header.
10 #include <top/modules/TOPLaserCalibrator/TOPLaserCalibratorModule.h>
11 
12 // TOP headers.
13 #include <top/modules/TOPLaserCalibrator/LaserCalibratorFit.h>
14 
15 // framework - DataStore
16 #include <framework/datastore/StoreArray.h>
17 
18 // framework aux
19 #include <framework/logging/Logger.h>
20 
21 // DataStore classes
22 //#include <framework/dataobjects/EventMetaData.h>
23 #include <top/dataobjects/TOPDigit.h>
24 
25 #include <TFile.h>
26 #include <TH1F.h>
27 #include <TMath.h>
28 #include <TTree.h>
29 #include <sstream>
30 
31 
32 using namespace std;
33 
34 namespace Belle2 {
39  //-----------------------------------------------------------------
41  //-----------------------------------------------------------------
42 
43  REG_MODULE(TOPLaserCalibrator);
44 
45  //-----------------------------------------------------------------
46  // Implementation
47  //-----------------------------------------------------------------
48 
49  TOPLaserCalibratorModule::TOPLaserCalibratorModule() : Module()
50  {
51  setDescription("T0 Calibration using the laser calibration system");
52 
53  std::vector<double> frange = {100, 0., 1.};
54  // Add parameters
55  addParam("dataFitOutput", m_dataFitOutput, "Output root file for data",
56  string("")); //output fitting results
57  addParam("mcInput", m_mcInput, "Input root file from MC", string(""));
58  addParam("channelT0constant", m_chT0C, "Output of channel T0 constant", string(""));
59  addParam("fitChannel", m_fitChannel, "set 0 - 511 to a specific channel; set 0 to all channels in the fit",
60  512);
61  addParam("barID", m_barID, "ID of TOP module to calibrate");
62  addParam("refCh", m_refCh, "reference channel of T0 constant");
63  addParam("fitMethod", m_fitMethod, "gauss: single gaussian; cb: single Crystal Ball; cb2: double Crystal Ball", string("gauss"));
64  addParam("fitRange", m_fitRange, "fit range[nbins, xmin, xmax]", frange);
65 
66  for (int i = 0; i < c_NumChannels; ++i) {
67  m_histo[i] = 0;
68  }
69  }
70 
72  {
73  }
74 
76  {
77  B2WARNING("You are using an old version of the laser fitter, now deprecated. This module has been superseded by the CAF collector TOPLaserCalibratorCollector and the CAF fitter TOPLocalCalFitter.");
78  m_digits.isRequired();
79  }
80 
82  {
83  }
84 
86  {
87  for (auto& digit : m_digits) {
88  if (m_barID != 0 && digit.getModuleID() != m_barID) continue; //if m_barID == 0, use all 16 slots(for MC test)
89  unsigned channel = digit.getChannel();
90  if (channel < c_NumChannels) {
91  auto histo = m_histo[channel];
92  if (!histo) {
93  stringstream ss;
94  ss << "chan" << channel ;
95  string name;
96  ss >> name;
97  string title = "Times " + name;
98  histo = new TH1F(name.c_str(), title.c_str(), (int)m_fitRange[0], m_fitRange[1], m_fitRange[2]);
99  m_histo[channel] = histo;
100  }
101  if (digit.getHitQuality() != TOPDigit::EHitQuality::c_Good) continue;
102  if (digit.getTime() < m_fitRange[1] || digit.getTime() > m_fitRange[2]) continue;
103  if (!digit.isTimeBaseCalibrated()) continue;
104  histo->Fill(digit.getTime()); //get Time from TOPDigit
105  }
106  }
107  }
108 
110  {
111  }
112 
114  {
115  std::vector<TH1F*> v_histo;
116  for (int i = 0; i < c_NumChannels; ++i) {
117  v_histo.push_back(m_histo[i]);
118  }
119 
120  double dataT[c_NumChannels] = {0};
121  double dataTErr[c_NumChannels] = {0};
122  double mcT[c_NumChannels] = {0};
123 
125  t0fit.setHist(v_histo); //set time hist vector
126  t0fit.setFitMethod(m_fitMethod);
127  t0fit.setFitRange(m_fitRange[1], m_fitRange[2]);
128  if (m_fitChannel == c_NumChannels) {
129  for (int i = 0; i < c_NumChannels; ++i) {
130  t0fit.fitChannel(i);
131  dataT[i] = t0fit.getFitT();
132  dataTErr[i] = t0fit.getFitTErr();
133  }
134  } else {
135  t0fit.fitChannel(m_fitChannel);
136  dataT[m_fitChannel] = t0fit.getFitT();
137  }
138  //cout<<"chisq=\t"<<t0fit.getFitChisq(6)<<endl;
139  t0fit.writeFile(m_dataFitOutput); //write to root file
140 
141  //read laser propagation time from MC (TOPChannelT0MC)
142  auto mcFile = new TFile(m_mcInput.c_str());
143  auto tree = static_cast<TTree*>(mcFile->Get("t0MC"));
144  int channel_mc;
145  double maxpos;
146 
147  tree->SetBranchAddress("channel", &channel_mc);
148  tree->SetBranchAddress("maxpos", &maxpos);
149  for (int i = 0; i < tree->GetEntries(); i++) {
150  tree->GetEntry(i);
151  mcT[channel_mc] = maxpos;
152  }
153 
154  //w.r.t. ref. channel
155  for (int i = 0; i < c_NumChannels; i++) {
156  if (i == m_refCh) continue;
157  dataT[i] = dataT[i] - dataT[m_refCh];
158  mcT[i] = mcT[i] - mcT[m_refCh];
159  }
160 
161  delete tree;
162  mcFile->Close();
163  delete mcFile;
164 
165  //calculate T0 Const, write to root file
166  auto file = new TFile(m_chT0C.c_str(), "RECREATE");
167  auto otree = new TTree("chT0", "Channel T0 Const");
168 
169  unsigned channel = 0;
170  double t0_const = 0;
171  double t0ConstRaw = 0;
172  double fittedTime = 0;
173  double fittedTimeError = 0;
174  double mcTime = 0;
175  double mcCorrection = 0;
176  double t0ConstError = 0;
177 
178  otree->Branch("fittedTime", &fittedTime, "fittedTime/D");
179  otree->Branch("fittedTimeError", &fittedTimeError, "fittedTimeError/D");
180  otree->Branch("mcTime", &mcTime, "mcTime/D");
181  otree->Branch("t0ConstRaw", &t0ConstRaw, "t0ConstRaw/D");
182  otree->Branch("mcCorrection", &mcCorrection, "mcCorrection/D");
183  otree->Branch("t0Const", &t0_const, "t0_const/D");
184  otree->Branch("t0ConstError", &t0ConstError, "t0ConstError/D");
185  otree->Branch("channel", &channel, "channel/I");
186  otree->Branch("slot", &m_barID, "slot/I");
187 
188  for (int i = 0; i < c_NumChannels; i++) {
189  if (i == m_refCh) {
190  fittedTime = dataT[i];
191  mcTime = mcT[i];
192  } else {
193  fittedTime = dataT[i] + dataT[m_refCh];
194  mcTime = mcT[i] + mcT[m_refCh];
195  }
196  fittedTimeError = dataTErr[i];
197  t0ConstRaw = dataT[i];
198  mcCorrection = mcT[i];
199  channel = i;
200  t0_const = dataT[i] - mcT[i];
201  // assuming that the MC error is negligible
202  t0ConstError = TMath::Sqrt(dataTErr[i] * dataTErr[i] + dataTErr[m_refCh] * dataTErr[m_refCh]);
203  if (i == m_refCh) {
204  t0_const = 0;
205  t0ConstError = dataTErr[m_refCh];
206  }
207  otree->Fill();
208  }
209  otree->Write();
210  delete otree;
211  file->Close();
212  delete file;
213 
214  // ...
215  // write to database; next to do ...
216  // ...
217  }
219 } // end Belle2 namespace
220 
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
std::vector< double > m_fitRange
fit range [nbins, xmin, xmax]
int m_fitChannel
set 0 - 511 to a specific pixelID in the fit; set 512 to fit all pixels in one slot
int m_barID
ID of TOP module to calibrate.
std::string m_mcInput
Input root file from MC.
StoreArray< TOPDigit > m_digits
collection of digits
int m_refCh
reference channel of T0 constant
std::string m_chT0C
Output of channel T0 constant.
TH1F * m_histo[c_NumChannels]
profile histograms
std::string m_fitMethod
gauss: single gaussian; cb: single Crystal Ball; cb2: double Crystal Ball
std::string m_dataFitOutput
output root file for data
A class do laser calibration fit provide different fitting method (under development)
void setHist(const std::vector< TH1F * > &hist)
set time hist of all channels in one moduleID
double getFitT()
get mean positon after fit
void setFitMethod(const std::string &method)
set time fit function
void setFitRange(double xmin=-200, double xmax=200)
set x range in the fit
int fitChannel(unsigned channel)
fit for a specific channel
double getFitTErr()
returns the error mean positon after fit
void writeFile(const std::string &outfile)
write fit result to a root file
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
Abstract base class for different kinds of events.