Belle II Software  release-08-00-10
eclWaveformTemplateCalibrationC2Algorithm.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 <ecl/calibration/eclWaveformTemplateCalibrationC2Algorithm.h>
11 
12 /* ECL headers. */
13 #include <ecl/dataobjects/ECLElementNumbers.h>
14 #include <ecl/digitization/OfflineFitFunction.h>
15 #include <ecl/dbobjects/ECLDigitWaveformParameters.h>
16 #include <ecl/dbobjects/ECLCrystalCalib.h>
17 
18 /* ROOT headers. */
19 #include <TFile.h>
20 #include <TGraph.h>
21 #include <TTree.h>
22 #include <TF1.h>
23 
24 #include<ctime>
25 
26 using namespace std;
27 using namespace Belle2;
28 using namespace ECL;
29 using namespace Calibration;
30 
32 eclWaveformTemplateCalibrationC2Algorithm::eclWaveformTemplateCalibrationC2Algorithm():
33  CalibrationAlgorithm("eclWaveformTemplateCalibrationC2Collector")
34 {
36  "Perform the photon template shape calibration using waveforms from high energy crystals from e+e- --> gamma gamma events"
37  );
38 
39 }
40 
41 
42 namespace {
43 
44  // Used to perform simultaneous fits of multiple waveforms
45  std::vector<TF1*> FitFunctions;
46  const double numberofADCPoints = 31.0;
47 
48  double fitf(double* x, double* par)
49  {
50 
51  double xtoeval = std::fmod(x[0], numberofADCPoints);
52  int whichFitFunctions = x[0] / numberofADCPoints;
53 
54  for (int i = 0; i < (int)FitFunctions.size(); i++) {
55  FitFunctions[i]->SetParameter(0, par[i]);
56  FitFunctions[i]->SetParameter(1, par[FitFunctions.size() + i]);
57  FitFunctions[i]->SetParameter(2, par[(2 * FitFunctions.size()) + i]);
58  FitFunctions[i]->FixParameter(3, 0);
59  for (int k = 0; k < 10; k++) {
60  FitFunctions[i]->SetParameter(4 + k, par[(3 * FitFunctions.size()) + k]);
61  }
62  }
63 
64  return FitFunctions[whichFitFunctions]->Eval(xtoeval * 0.5);
65  }
66 
67 }
68 
70 {
71 
72  B2INFO("Reading ECLCrystalCalib payload: eclWaveformTemplateCalibrationC1MaxResLimit");
73  DBObjPtr<ECLCrystalCalib> existingeclWaveformTemplateCalibrationC1MaxResLimit("eclWaveformTemplateCalibrationC1MaxResLimit");
74  auto runs = getRunList();
75  ExpRun chosenRun = runs.front();
76  // After here your DBObjPtrs are correct
77  updateDBObjPtrs(1, chosenRun.second, chosenRun.first);
78 
80  gROOT->SetBatch();
81 
82  std::vector<double> cellIDArray;
83  std::vector<double> maxResidualArray; // used to quantify fit result
84  std::vector<double> limitResidualArray; // what was final resLimit used
85  std::vector<double> parLimitFactorArray; // what was final parLimitFactorArray used
86 
88  TFile* histfile = new TFile(m_outputName.c_str(), "recreate");
89 
91  TFile* f_PhotonTemplateOutput = new TFile(Form("PhotonShapes_Low%d_High%d.root", m_firstCellID, m_lastCellID), "RECREATE");
92  TTree* mtree = new TTree("mtree", "");
93  std::vector<double> PhotonWaveformArray(100000);
94  mtree->Branch("PhotonArray", PhotonWaveformArray.data(), "PhotonWaveformArray[100000]/D");
95 
97  auto tree = getObjectPtr<TTree>("tree");
98  int CellID;
99  tree->SetBranchAddress("CellID", &CellID);
100  std::vector<int> Waveform(m_NumberofADCPoints);
101  std::vector<int> XValues(m_NumberofADCPoints);
102  for (int i = 0; i < m_NumberofADCPoints; i++) {
103  tree->SetBranchAddress(Form("ADC%d", i), &Waveform[i]);
104  XValues[i] = i;
105  }
106 
107  std::time_t t = std::time(0);
108 
110  int AttemptCounter = 0;
111 
113  ECLDigitWaveformParameters* PhotonParameters = new ECLDigitWaveformParameters();
114 
117  double ParMin11t[11];
118 
120  for (int CellID_i = m_firstCellID; CellID_i <= m_lastCellID; CellID_i++) {
121 
123  if (CellID_i > 7776 || CellID_i < 1153) {
124  ParMin11t[0] = 20.3216;
125  ParMin11t[1] = -0.0206266;
126  ParMin11t[2] = 0.313928;
127  ParMin11t[3] = 0.589646;
128  ParMin11t[4] = 0.455526;
129  ParMin11t[5] = 1.03656;
130  ParMin11t[6] = 0.000822467;
131  ParMin11t[7] = 45.1574;
132  ParMin11t[8] = 0.716034;
133  ParMin11t[9] = 0.616753;
134  ParMin11t[10] = 0.0851222;
135  } else {
136  ParMin11t[0] = 24.6176;
137  ParMin11t[1] = 0.00725002;
138  ParMin11t[2] = 0.601578;
139  ParMin11t[3] = 0.491976;
140  ParMin11t[4] = 0.601034;
141  ParMin11t[5] = 0.601684;
142  ParMin11t[6] = -0.0103788;
143  ParMin11t[7] = 2.22615;
144  ParMin11t[8] = 0.671294;
145  ParMin11t[9] = 0.529878;
146  ParMin11t[10] = 0.0757927;
147  }
148 
149  double resLimit = 2 * existingeclWaveformTemplateCalibrationC1MaxResLimit->getCalibVector()[CellID_i -
150  1];
151  double resLimitOriginal = resLimit;
152 
154  std::vector<int> EntriesToSkip;
155 
157  double maxResidual = 1000.0;
158 
160  bool PASS = false;
161  while (PASS == false) {
162 
164  std::vector<double> xValuesToFit;
165  std::vector<double> yValuesToFit;
166 
168  std::vector<double> guessBaseline;
169  std::vector<double> guessAmp;
170  std::vector<double> guessTime;
171 
173  std::vector<int> NtupleEntries;
174 
175  int counter = 0; // counts entry number in xValuesToFit
176  int counterWaveforms = 0; // counts numver of waveforms selected
177 
178  for (int i = 0; i < tree->GetEntries(); i++) {
179 
181  bool skipEvent = false;
182  for (int k = 0; k < (int)EntriesToSkip.size(); k++) {
183  if (EntriesToSkip[k] == i) skipEvent = true;
184  }
185  if (skipEvent) continue;
186 
187  tree->GetEntry(i);
188 
189  if (CellID != CellID_i) continue;
190 
191  double maxval = 0;
192  double maxIndex = 0;
193  for (int j = 0; j < m_NumberofADCPoints; j++) {
194  xValuesToFit.push_back(counter);
195  yValuesToFit.push_back(Waveform[j]);
196  if (Waveform[j] > maxval) {
197  maxval = Waveform[j];
198  maxIndex = j;
199  }
200  counter++;
201  }
202 
204  guessBaseline.push_back(Waveform[0]);
205  guessAmp.push_back(maxval);
206  guessTime.push_back((maxIndex - 4.5) * 0.5);
208  NtupleEntries.push_back(i);
209  B2INFO("Entry: " << i);
210 
211  counterWaveforms++;
212 
213  if (counterWaveforms == m_CollectorLimit) break;
214 
215  }
216 
218  B2INFO("CellID " << CellID_i << " counterWaveforms = " << counterWaveforms);
219 
220  if (counterWaveforms < m_TotalCountsThreshold) {
222  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: warning total entries for cell ID " << CellID_i << " is only: " <<
223  counterWaveforms << " Requirement is : " << m_TotalCountsThreshold);
224  EntriesToSkip.clear();
225  resLimit *= 2;
226  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: warning " << CellID_i << " resLimit is doubled to" << resLimit << " start was "
227  << resLimitOriginal);
228  }
229 
230 
232  auto gWaveformToFit = new TGraph(xValuesToFit.size(), xValuesToFit.data(), yValuesToFit.data());
233  gWaveformToFit->SetName(Form("gWaveformToFit_%d", int(CellID_i)));
234 
238  FitFunctions.clear();
239  for (int i = 0; i < counterWaveforms; i++) {
240  FitFunctions.push_back(new TF1(Form("Shp_%d", i), Belle2::ECL::WaveFuncTwoComponent, 0, 30.5, 26));
241  FitFunctions[i]->SetNpx(10000);
242  FitFunctions[i]->FixParameter(3, 0);
243  for (int k = 0; k < 10; k++) {
244  FitFunctions[i]->SetParameter(4 + k, ParMin11t[k + 1]);
245  FitFunctions[i]->FixParameter(10 + 4 + k, ParMin11t[k + 1]);
246  }
247  FitFunctions[i]->FixParameter(24, ParMin11t[0]);
248  FitFunctions[i]->FixParameter(25, 1);
249  }
250 
252  TF1* TotalFitFunction = new TF1("TotalFitFunction", fitf, 0, counterWaveforms * m_NumberofADCPoints,
253  (3 * FitFunctions.size()) + 10);
254 
256  int FFsize = FitFunctions.size();
257  for (int i = 0; i < FFsize; i++) {
258  TotalFitFunction->SetParameter(i, guessTime[i]);
259  TotalFitFunction->SetParameter(FFsize + i, guessBaseline[i]);
260  TotalFitFunction->SetParameter((2 * FFsize) + i, guessAmp[i]);
261  for (int k = 0; k < 10; k++) {
262  TotalFitFunction->SetParameter((3 * FFsize) + k, ParMin11t[k + 1]);
263  if (m_ParamLimitFactor < 2) {
264  TotalFitFunction->SetParLimits((3 * FFsize) + k, ParMin11t[k + 1] - m_ParamLimitFactor * fabs(ParMin11t[k + 1]),
265  ParMin11t[k + 1] + m_ParamLimitFactor * fabs(ParMin11t[k + 1]));
266  } else {
267  TotalFitFunction->ReleaseParameter((3 * FFsize) + k);
268  }
269  }
270  }
271 
273  gWaveformToFit->Fit("TotalFitFunction", "Q M W N 0 R", "", 0, counterWaveforms * m_NumberofADCPoints);
274 
276  std::vector<int> FitResultY;
277  std::vector<int> FitResultX;
278  int maxResidualWaveformID = 0; // Used to remove waveforms with potential pile-up outside baseline
279 
281  maxResidual = 0.0;
282  double npts = xValuesToFit.size();
283  double maxResidualOld = 0.0;
284  for (int k = 0; k < npts; k++) {
285  double xVal = xValuesToFit[k];
286  double yVal = TotalFitFunction->Eval(xVal);
287  FitResultX.push_back(xVal);
288  FitResultY.push_back(yVal);
289  double diff = fabs(yValuesToFit[k] - yVal);
290  if (diff > maxResidual) {
291  maxResidual = diff;
292  maxResidualWaveformID = (k / m_NumberofADCPoints);
293  maxResidualOld = fabs(yValuesToFit[k] / yVal);
294  }
295  }
296 
297  // Checking if fit matches the data.
298  if (maxResidual > resLimit) {
299 
300  B2INFO("FAIL: CellID_i " << CellID_i << " maxResidual " << maxResidual << " removing entry: " <<
301  NtupleEntries[maxResidualWaveformID] <<
302  " which was waveform number " << maxResidualWaveformID << " resLimit was " << resLimit << " , resLimit started at " <<
303  resLimitOriginal);
304  B2INFO("Old maxResidual of Data/Fit was " << maxResidualOld);
305 
306  B2INFO("Iter Time = " << std::time(0) - t << std::endl);
307  t = std::time(0);
308 
309  std::cout << "FAIL: CellID_i " << CellID_i << " maxResidual " << maxResidual << " removing entry: " <<
310  NtupleEntries[maxResidualWaveformID] <<
311  " which was waveform number " << maxResidualWaveformID << " resLimit was " << resLimit << " , resLimit started at " <<
312  resLimitOriginal << std::endl;
313  std::cout << "wave = [";
314  for (int k = 0; k < npts; k++) {
315  std::cout << yValuesToFit[k];
316  if (k < (npts - 1)) {
317  std::cout << ",";
318  } else {
319  std::cout << "]" << std::endl;
320  }
321  }
322  std::cout << "fitRes = [";
323  for (int k = 0; k < npts; k++) {
324  std::cout << TotalFitFunction->Eval(xValuesToFit[k]);
325  if (k < (npts - 1)) {
326  std::cout << ",";
327  } else {
328  std::cout << "]" << std::endl;
329  }
330  }
331 
333  EntriesToSkip.push_back(NtupleEntries[maxResidualWaveformID]);
334 
335  AttemptCounter++;
336 
338  if (counterWaveforms < m_SimutaniousFitLimit) AttemptCounter = m_AttemptLimit;
339 
341  if (AttemptCounter == m_AttemptLimit) {
342 
344 
345  B2INFO("AttemptCounter reach limit: " << AttemptCounter << " counterWaveforms: " << counterWaveforms);
346  B2INFO("Increasing m_ParamLimitFactor to " << m_ParamLimitFactor);
347 
349  EntriesToSkip.clear();
350  AttemptCounter = 0;
351 
354  resLimit *= m_ResLimitIterator;
355  B2INFO("Increasing resLimit to " << resLimit);
357  }
358  }
359 
360  } else {
361 
362  B2INFO("PASS: CellID_i " << CellID_i << " maxResidual " << maxResidual << " number of waveforms used was " << counterWaveforms <<
363  " resLimit was " << resLimit);
364 
365  PASS = true;
366 
367  limitResidualArray.push_back(resLimit);
368  parLimitFactorArray.push_back(m_ParamLimitFactor);
369 
371  AttemptCounter = 0;
373 
374  auto gFitResult = new TGraph(FitResultX.size(), FitResultX.data(), FitResultY.data());
375  gFitResult->SetName(Form("gFitResult_%d", int(CellID_i)));
376 
378  cellIDArray.push_back(CellID_i);
379  maxResidualArray.push_back(maxResidual);
380 
382  histfile->cd();
383  gWaveformToFit->Write();
384  gFitResult->Write();
385 
387  float tempPhotonPar11[11];
388  tempPhotonPar11[0] = ParMin11t[0];
389  for (unsigned int k = 0; k < 10; k++) tempPhotonPar11[k + 1] = TotalFitFunction->GetParameter((3 * FFsize) + k);
390 
392  FitFunctions[0]->SetParameter(0, 0);
393  FitFunctions[0]->SetParameter(1, 0);
394  FitFunctions[0]->SetParameter(2, 1);
395  for (int k = 0; k < 10; k++) {
396  FitFunctions[0]->SetParameter(4 + k, tempPhotonPar11[k + 1]);
397  FitFunctions[0]->SetParameter(10 + 4 + k, tempPhotonPar11[k + 1]);
398  }
399  FitFunctions[0]->FixParameter(24, ParMin11t[0]);
400  FitFunctions[0]->FixParameter(25, 1);
401 
403  double MaxVal = -1.0;
404  const double cnpts = 2000;
405  for (int k = 0; k < cnpts; k++) {
406  double xVal = (k * double(m_NumberofADCPoints) / cnpts);
407  double yVal = FitFunctions[0]->Eval(xVal);
408  if (yVal > MaxVal) MaxVal = yVal;
409  }
410  B2INFO("MaxVal " << MaxVal);
411  tempPhotonPar11[0] /= MaxVal;
412  FitFunctions[0]->FixParameter(24, tempPhotonPar11[0]);
413 
415  PhotonParameters->setTemplateParameters(CellID_i, tempPhotonPar11, tempPhotonPar11, tempPhotonPar11);
416 
418  for (unsigned int k = 0; k < PhotonWaveformArray.size();
419  k++) PhotonWaveformArray[k] = FitFunctions[0]->Eval(((double)k) * (1. / 1000.)) ;
420  mtree->Fill();
421 
422  }
423  for (int w = 0; w < (int)FitFunctions.size(); w++) FitFunctions[w]->Delete();
424  TotalFitFunction->Delete() ;
425  gWaveformToFit->Delete();
426  }
427  }
428 
430  histfile->cd();
431  auto gmaxResidual = new TGraph(cellIDArray.size(), cellIDArray.data(), maxResidualArray.data());
432  gmaxResidual->SetName("gmaxResidual");
433  auto glimitResidualArray = new TGraph(cellIDArray.size(), cellIDArray.data(), limitResidualArray.data());
434  glimitResidualArray->SetName("glimitResidualArray");
435  auto gparLimitFactorArray = new TGraph(cellIDArray.size(), cellIDArray.data(), parLimitFactorArray.data());
436  gparLimitFactorArray->SetName("gparLimitFactorArray");
437 
438  gmaxResidual->Write();
439  glimitResidualArray->Write();
440  gparLimitFactorArray->Write();
441  histfile->Write();
442  histfile->Close();
443  delete histfile;
444 
445  f_PhotonTemplateOutput->cd();
446  mtree->Write();
447  f_PhotonTemplateOutput->Write();
448  f_PhotonTemplateOutput->Close();
449  delete f_PhotonTemplateOutput;
450 
452  saveCalibration(PhotonParameters, Form("PhotonParameters_CellID%d_CellID%d", m_firstCellID, m_lastCellID));
453  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: successfully stored " << Form("PhotonParameters_CellID%d_CellID%d",
454  m_firstCellID, m_lastCellID) << " constants");
455 
456  return c_OK;
457 }
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
DB object to store photon, hadron and diode shape parameters.
void setTemplateParameters(int cellID, const float photonInput[11], const float hadronInput[11], const float diodeInput[11])
Set photon, hadron and diode template parameters for crystal.
int m_AttemptLimit
Number of attempts before increasing parameter limits or resLimt.
double m_ParamLimitFactor
Factor to determine parameter limits in fit.
const int m_SimutaniousFitLimit
Min number waveforms required for simultaneous fit.
const int m_TotalCountsThreshold
Min number waveforms required per crystal.
Abstract base class for different kinds of events.
Struct containing exp number and run number.
Definition: Splitter.h:51