Belle II Software  release-08-01-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 /* Basf2 headers. */
19 #include <framework/database/DBObjPtr.h>
20 
21 /* ROOT headers. */
22 #include <TFile.h>
23 #include <TGraph.h>
24 #include <TTree.h>
25 #include <TF1.h>
26 
27 /* C++ headers. */
28 #include <ctime>
29 
30 using namespace std;
31 using namespace Belle2;
32 using namespace ECL;
33 using namespace Calibration;
34 
36 eclWaveformTemplateCalibrationC2Algorithm::eclWaveformTemplateCalibrationC2Algorithm():
37  CalibrationAlgorithm("eclWaveformTemplateCalibrationC2Collector")
38 {
40  "Perform the photon template shape calibration using waveforms from high energy crystals from e+e- --> gamma gamma events"
41  );
42 
43 }
44 
45 
46 namespace {
47 
48  // Used to perform simultaneous fits of multiple waveforms
49  std::vector<TF1*> FitFunctions;
50  const double numberofADCPoints = 31.0;
51 
52  double fitf(double* x, double* par)
53  {
54 
55  double xtoeval = std::fmod(x[0], numberofADCPoints);
56  int whichFitFunctions = x[0] / numberofADCPoints;
57 
58  for (int i = 0; i < (int)FitFunctions.size(); i++) {
59  FitFunctions[i]->SetParameter(0, par[i]);
60  FitFunctions[i]->SetParameter(1, par[FitFunctions.size() + i]);
61  FitFunctions[i]->SetParameter(2, par[(2 * FitFunctions.size()) + i]);
62  FitFunctions[i]->FixParameter(3, 0);
63  for (int k = 0; k < 10; k++) {
64  FitFunctions[i]->SetParameter(4 + k, par[(3 * FitFunctions.size()) + k]);
65  }
66  }
67 
68  return FitFunctions[whichFitFunctions]->Eval(xtoeval * 0.5);
69  }
70 
71 }
72 
74 {
75 
76  B2INFO("Reading ECLCrystalCalib payload: eclWaveformTemplateCalibrationC1MaxResLimit");
77  DBObjPtr<ECLCrystalCalib> existingeclWaveformTemplateCalibrationC1MaxResLimit("eclWaveformTemplateCalibrationC1MaxResLimit");
78  auto runs = getRunList();
79  ExpRun chosenRun = runs.front();
80  // After here your DBObjPtrs are correct
81  updateDBObjPtrs(1, chosenRun.second, chosenRun.first);
82 
84  gROOT->SetBatch();
85 
86  std::vector<double> cellIDArray;
87  std::vector<double> maxResidualArray; // used to quantify fit result
88  std::vector<double> limitResidualArray; // what was final resLimit used
89  std::vector<double> parLimitFactorArray; // what was final parLimitFactorArray used
90 
92  TFile* histfile = new TFile(m_outputName.c_str(), "recreate");
93 
95  TFile* f_PhotonTemplateOutput = new TFile(Form("PhotonShapes_Low%d_High%d.root", m_firstCellID, m_lastCellID), "RECREATE");
96  TTree* mtree = new TTree("mtree", "");
97  std::vector<double> PhotonWaveformArray(100000);
98  mtree->Branch("PhotonArray", PhotonWaveformArray.data(), "PhotonWaveformArray[100000]/D");
99 
101  auto tree = getObjectPtr<TTree>("tree");
102  int CellID;
103  tree->SetBranchAddress("CellID", &CellID);
104  std::vector<int> Waveform(m_NumberofADCPoints);
105  std::vector<int> XValues(m_NumberofADCPoints);
106  for (int i = 0; i < m_NumberofADCPoints; i++) {
107  tree->SetBranchAddress(Form("ADC%d", i), &Waveform[i]);
108  XValues[i] = i;
109  }
110 
111  std::time_t t = std::time(0);
112 
114  int AttemptCounter = 0;
115 
117  ECLDigitWaveformParameters* PhotonParameters = new ECLDigitWaveformParameters();
118 
121  double ParMin11t[11];
122 
124  for (int CellID_i = m_firstCellID; CellID_i <= m_lastCellID; CellID_i++) {
125 
127  if (CellID_i > 7776 || CellID_i < 1153) {
128  ParMin11t[0] = 20.3216;
129  ParMin11t[1] = -0.0206266;
130  ParMin11t[2] = 0.313928;
131  ParMin11t[3] = 0.589646;
132  ParMin11t[4] = 0.455526;
133  ParMin11t[5] = 1.03656;
134  ParMin11t[6] = 0.000822467;
135  ParMin11t[7] = 45.1574;
136  ParMin11t[8] = 0.716034;
137  ParMin11t[9] = 0.616753;
138  ParMin11t[10] = 0.0851222;
139  } else {
140  ParMin11t[0] = 24.6176;
141  ParMin11t[1] = 0.00725002;
142  ParMin11t[2] = 0.601578;
143  ParMin11t[3] = 0.491976;
144  ParMin11t[4] = 0.601034;
145  ParMin11t[5] = 0.601684;
146  ParMin11t[6] = -0.0103788;
147  ParMin11t[7] = 2.22615;
148  ParMin11t[8] = 0.671294;
149  ParMin11t[9] = 0.529878;
150  ParMin11t[10] = 0.0757927;
151  }
152 
153  double resLimit = 2 * existingeclWaveformTemplateCalibrationC1MaxResLimit->getCalibVector()[CellID_i -
154  1];
155  double resLimitOriginal = resLimit;
156 
158  std::vector<int> EntriesToSkip;
159 
161  double maxResidual = 1000.0;
162 
164  bool PASS = false;
165  while (PASS == false) {
166 
168  std::vector<double> xValuesToFit;
169  std::vector<double> yValuesToFit;
170 
172  std::vector<double> guessBaseline;
173  std::vector<double> guessAmp;
174  std::vector<double> guessTime;
175 
177  std::vector<int> NtupleEntries;
178 
179  int counter = 0; // counts entry number in xValuesToFit
180  int counterWaveforms = 0; // counts numver of waveforms selected
181 
182  for (int i = 0; i < tree->GetEntries(); i++) {
183 
185  bool skipEvent = false;
186  for (int k = 0; k < (int)EntriesToSkip.size(); k++) {
187  if (EntriesToSkip[k] == i) skipEvent = true;
188  }
189  if (skipEvent) continue;
190 
191  tree->GetEntry(i);
192 
193  if (CellID != CellID_i) continue;
194 
195  double maxval = 0;
196  double maxIndex = 0;
197  for (int j = 0; j < m_NumberofADCPoints; j++) {
198  xValuesToFit.push_back(counter);
199  yValuesToFit.push_back(Waveform[j]);
200  if (Waveform[j] > maxval) {
201  maxval = Waveform[j];
202  maxIndex = j;
203  }
204  counter++;
205  }
206 
208  guessBaseline.push_back(Waveform[0]);
209  guessAmp.push_back(maxval);
210  guessTime.push_back((maxIndex - 4.5) * 0.5);
212  NtupleEntries.push_back(i);
213  B2INFO("Entry: " << i);
214 
215  counterWaveforms++;
216 
217  if (counterWaveforms == m_CollectorLimit) break;
218 
219  }
220 
222  B2INFO("CellID " << CellID_i << " counterWaveforms = " << counterWaveforms);
223 
224  if (counterWaveforms < m_TotalCountsThreshold) {
226  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: warning total entries for cell ID " << CellID_i << " is only: " <<
227  counterWaveforms << " Requirement is : " << m_TotalCountsThreshold);
228  EntriesToSkip.clear();
229  resLimit *= 2;
230  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: warning " << CellID_i << " resLimit is doubled to" << resLimit << " start was "
231  << resLimitOriginal);
232  }
233 
234 
236  auto gWaveformToFit = new TGraph(xValuesToFit.size(), xValuesToFit.data(), yValuesToFit.data());
237  gWaveformToFit->SetName(Form("gWaveformToFit_%d", int(CellID_i)));
238 
242  FitFunctions.clear();
243  for (int i = 0; i < counterWaveforms; i++) {
244  FitFunctions.push_back(new TF1(Form("Shp_%d", i), Belle2::ECL::WaveFuncTwoComponent, 0, 30.5, 26));
245  FitFunctions[i]->SetNpx(10000);
246  FitFunctions[i]->FixParameter(3, 0);
247  for (int k = 0; k < 10; k++) {
248  FitFunctions[i]->SetParameter(4 + k, ParMin11t[k + 1]);
249  FitFunctions[i]->FixParameter(10 + 4 + k, ParMin11t[k + 1]);
250  }
251  FitFunctions[i]->FixParameter(24, ParMin11t[0]);
252  FitFunctions[i]->FixParameter(25, 1);
253  }
254 
256  TF1* TotalFitFunction = new TF1("TotalFitFunction", fitf, 0, counterWaveforms * m_NumberofADCPoints,
257  (3 * FitFunctions.size()) + 10);
258 
260  int FFsize = FitFunctions.size();
261  for (int i = 0; i < FFsize; i++) {
262  TotalFitFunction->SetParameter(i, guessTime[i]);
263  TotalFitFunction->SetParameter(FFsize + i, guessBaseline[i]);
264  TotalFitFunction->SetParameter((2 * FFsize) + i, guessAmp[i]);
265  for (int k = 0; k < 10; k++) {
266  TotalFitFunction->SetParameter((3 * FFsize) + k, ParMin11t[k + 1]);
267  if (m_ParamLimitFactor < 2) {
268  TotalFitFunction->SetParLimits((3 * FFsize) + k, ParMin11t[k + 1] - m_ParamLimitFactor * fabs(ParMin11t[k + 1]),
269  ParMin11t[k + 1] + m_ParamLimitFactor * fabs(ParMin11t[k + 1]));
270  } else {
271  TotalFitFunction->ReleaseParameter((3 * FFsize) + k);
272  }
273  }
274  }
275 
277  gWaveformToFit->Fit("TotalFitFunction", "Q M W N 0 R", "", 0, counterWaveforms * m_NumberofADCPoints);
278 
280  std::vector<int> FitResultY;
281  std::vector<int> FitResultX;
282  int maxResidualWaveformID = 0; // Used to remove waveforms with potential pile-up outside baseline
283 
285  maxResidual = 0.0;
286  double npts = xValuesToFit.size();
287  double maxResidualOld = 0.0;
288  for (int k = 0; k < npts; k++) {
289  double xVal = xValuesToFit[k];
290  double yVal = TotalFitFunction->Eval(xVal);
291  FitResultX.push_back(xVal);
292  FitResultY.push_back(yVal);
293  double diff = fabs(yValuesToFit[k] - yVal);
294  if (diff > maxResidual) {
295  maxResidual = diff;
296  maxResidualWaveformID = (k / m_NumberofADCPoints);
297  maxResidualOld = fabs(yValuesToFit[k] / yVal);
298  }
299  }
300 
301  // Checking if fit matches the data.
302  if (maxResidual > resLimit) {
303 
304  B2INFO("FAIL: CellID_i " << CellID_i << " maxResidual " << maxResidual << " removing entry: " <<
305  NtupleEntries[maxResidualWaveformID] <<
306  " which was waveform number " << maxResidualWaveformID << " resLimit was " << resLimit << " , resLimit started at " <<
307  resLimitOriginal);
308  B2INFO("Old maxResidual of Data/Fit was " << maxResidualOld);
309 
310  B2INFO("Iter Time = " << std::time(0) - t << std::endl);
311  t = std::time(0);
312 
313  std::cout << "FAIL: CellID_i " << CellID_i << " maxResidual " << maxResidual << " removing entry: " <<
314  NtupleEntries[maxResidualWaveformID] <<
315  " which was waveform number " << maxResidualWaveformID << " resLimit was " << resLimit << " , resLimit started at " <<
316  resLimitOriginal << std::endl;
317  std::cout << "wave = [";
318  for (int k = 0; k < npts; k++) {
319  std::cout << yValuesToFit[k];
320  if (k < (npts - 1)) {
321  std::cout << ",";
322  } else {
323  std::cout << "]" << std::endl;
324  }
325  }
326  std::cout << "fitRes = [";
327  for (int k = 0; k < npts; k++) {
328  std::cout << TotalFitFunction->Eval(xValuesToFit[k]);
329  if (k < (npts - 1)) {
330  std::cout << ",";
331  } else {
332  std::cout << "]" << std::endl;
333  }
334  }
335 
337  EntriesToSkip.push_back(NtupleEntries[maxResidualWaveformID]);
338 
339  AttemptCounter++;
340 
342  if (counterWaveforms < m_SimutaniousFitLimit) AttemptCounter = m_AttemptLimit;
343 
345  if (AttemptCounter == m_AttemptLimit) {
346 
348 
349  B2INFO("AttemptCounter reach limit: " << AttemptCounter << " counterWaveforms: " << counterWaveforms);
350  B2INFO("Increasing m_ParamLimitFactor to " << m_ParamLimitFactor);
351 
353  EntriesToSkip.clear();
354  AttemptCounter = 0;
355 
358  resLimit *= m_ResLimitIterator;
359  B2INFO("Increasing resLimit to " << resLimit);
361  }
362  }
363 
364  } else {
365 
366  B2INFO("PASS: CellID_i " << CellID_i << " maxResidual " << maxResidual << " number of waveforms used was " << counterWaveforms <<
367  " resLimit was " << resLimit);
368 
369  PASS = true;
370 
371  limitResidualArray.push_back(resLimit);
372  parLimitFactorArray.push_back(m_ParamLimitFactor);
373 
375  AttemptCounter = 0;
377 
378  auto gFitResult = new TGraph(FitResultX.size(), FitResultX.data(), FitResultY.data());
379  gFitResult->SetName(Form("gFitResult_%d", int(CellID_i)));
380 
382  cellIDArray.push_back(CellID_i);
383  maxResidualArray.push_back(maxResidual);
384 
386  histfile->cd();
387  gWaveformToFit->Write();
388  gFitResult->Write();
389 
391  float tempPhotonPar11[11];
392  tempPhotonPar11[0] = ParMin11t[0];
393  for (unsigned int k = 0; k < 10; k++) tempPhotonPar11[k + 1] = TotalFitFunction->GetParameter((3 * FFsize) + k);
394 
396  FitFunctions[0]->SetParameter(0, 0);
397  FitFunctions[0]->SetParameter(1, 0);
398  FitFunctions[0]->SetParameter(2, 1);
399  for (int k = 0; k < 10; k++) {
400  FitFunctions[0]->SetParameter(4 + k, tempPhotonPar11[k + 1]);
401  FitFunctions[0]->SetParameter(10 + 4 + k, tempPhotonPar11[k + 1]);
402  }
403  FitFunctions[0]->FixParameter(24, ParMin11t[0]);
404  FitFunctions[0]->FixParameter(25, 1);
405 
407  double MaxVal = -1.0;
408  const double cnpts = 2000;
409  for (int k = 0; k < cnpts; k++) {
410  double xVal = (k * double(m_NumberofADCPoints) / cnpts);
411  double yVal = FitFunctions[0]->Eval(xVal);
412  if (yVal > MaxVal) MaxVal = yVal;
413  }
414  B2INFO("MaxVal " << MaxVal);
415  tempPhotonPar11[0] /= MaxVal;
416  FitFunctions[0]->FixParameter(24, tempPhotonPar11[0]);
417 
419  PhotonParameters->setTemplateParameters(CellID_i, tempPhotonPar11, tempPhotonPar11, tempPhotonPar11);
420 
422  for (unsigned int k = 0; k < PhotonWaveformArray.size();
423  k++) PhotonWaveformArray[k] = FitFunctions[0]->Eval(((double)k) * (1. / 1000.)) ;
424  mtree->Fill();
425 
426  }
427  for (int w = 0; w < (int)FitFunctions.size(); w++) FitFunctions[w]->Delete();
428  TotalFitFunction->Delete() ;
429  gWaveformToFit->Delete();
430  }
431  }
432 
434  histfile->cd();
435  auto gmaxResidual = new TGraph(cellIDArray.size(), cellIDArray.data(), maxResidualArray.data());
436  gmaxResidual->SetName("gmaxResidual");
437  auto glimitResidualArray = new TGraph(cellIDArray.size(), cellIDArray.data(), limitResidualArray.data());
438  glimitResidualArray->SetName("glimitResidualArray");
439  auto gparLimitFactorArray = new TGraph(cellIDArray.size(), cellIDArray.data(), parLimitFactorArray.data());
440  gparLimitFactorArray->SetName("gparLimitFactorArray");
441 
442  gmaxResidual->Write();
443  glimitResidualArray->Write();
444  gparLimitFactorArray->Write();
445  histfile->Write();
446  histfile->Close();
447  delete histfile;
448 
449  f_PhotonTemplateOutput->cd();
450  mtree->Write();
451  f_PhotonTemplateOutput->Write();
452  f_PhotonTemplateOutput->Close();
453  delete f_PhotonTemplateOutput;
454 
456  saveCalibration(PhotonParameters, Form("PhotonParameters_CellID%d_CellID%d", m_firstCellID, m_lastCellID));
457  B2INFO("eclWaveformTemplateCalibrationC2Algorithm: successfully stored " << Form("PhotonParameters_CellID%d_CellID%d",
458  m_firstCellID, m_lastCellID) << " constants");
459 
460  return c_OK;
461 }
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