Belle II Software  release-06-02-00
eclMergingCrystalTimingAlgorithm.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 
10 #include <ecl/calibration/eclMergingCrystalTimingAlgorithm.h>
11 #include <ecl/dbobjects/ECLCrystalCalib.h>
12 #include <ecl/dbobjects/ECLReferenceCrystalPerCrateCalib.h>
13 #include <ecl/utility/ECLChannelMapper.h>
14 #include <ecl/digitization/EclConfiguration.h>
15 #include "TH1F.h"
16 #include "TString.h"
17 #include "TFile.h"
18 #include "TDirectory.h"
19 
20 using namespace std;
21 using namespace Belle2;
22 using namespace ECL;
23 using namespace Calibration;
24 
26 eclMergingCrystalTimingAlgorithm::eclMergingCrystalTimingAlgorithm():
27  CalibrationAlgorithm("DummyCollector"),
28  readPrevCrysPayload(false),
29  m_ECLCrystalTimeOffsetBhabha("ECLCrystalTimeOffsetBhabha"),
30  m_ECLCrystalTimeOffsetBhabhaGamma("ECLCrystalTimeOffsetBhabhaGamma"),
31  m_ECLCrystalTimeOffsetCosmic("ECLCrystalTimeOffsetCosmic"),
32  m_ECLReferenceCrystalPerCrateCalib("ECLReferenceCrystalPerCrateCalib"),
33  m_ECLCrateTimeOffset("ECLCrateTimeOffset")
34 {
36  "Perform time calibration of ecl crystals by combining previous values from the DB for different calibrations."
37  );
38 }
39 
41 {
43  B2INFO("readPrevCrysPayload = " << readPrevCrysPayload);
44 
45  //------------------------------------------------------------------------
46  // Get the input run list to use to update the DBObjectPtrs
47  auto runs = getRunList();
48  /* Take the first run. For the crystal cosmic calibrations, because of the crate
49  calibrations, there is not a known correct run to use within the range. */
50  ExpRun chosenRun = runs.front();
51  B2INFO("merging using the ExpRun (" << chosenRun.second << "," << chosenRun.first << ")");
52  // After here your DBObjPtrs are correct
53  updateDBObjPtrs(1, chosenRun.second, chosenRun.first);
54 
55 
56  //------------------------------------------------------------------------
57  /* With the database set up, load the ecl channel mapper to access the
58  crystal to crate mapping information */
59  unique_ptr<ECLChannelMapper> crystalMapper(new ECL::ECLChannelMapper());
60  crystalMapper->initFromDB();
61 
62  //------------------------------------------------------------------------
63  // Test the DBObjects we want to exist and fail if not all of them do.
64  bool allObjectsFound = true;
65 
67  // Check that the payloads we want to merge are sufficiently loaded
69  allObjectsFound = false;
70  B2ERROR("No valid DBObject found for 'ECLCrystalTimeOffsetBhabha'");
71  }
73  allObjectsFound = false;
74  B2ERROR("No valid DBObject found for 'ECLCrystalTimeOffsetBhabhaGamma'");
75  }
77  allObjectsFound = false;
78  B2ERROR("No valid DBObject found for 'ECLCrystalTimeOffsetCosmic'");
79  }
80 
81  // Check that the reference crystal payload is loaded
83  allObjectsFound = false;
84  B2ERROR("No valid DBObject found for 'ECLReferenceCrystalPerCrateCalib'");
85  }
86 
87  // Check that the crate payload is loaded (used for transforming cosmic payload)
88  if (!m_ECLCrateTimeOffset) {
89  allObjectsFound = false;
90  B2ERROR("No valid DBObject found for 'ECLCrateTimeOffset'");
91  }
92 
93 
94  if (allObjectsFound) {
95  B2INFO("Valid objects found for 'ECLCrystalTimeOffsetBhabha', 'ECLCrystalTimeOffsetBhabhaGamma' and 'ECLCrystalTimeOffsetCosmic')");
96  B2INFO("Valid object found for 'ECLReferenceCrystalPerCrateCalib'");
97  B2INFO("Valid object found for 'ECLCrateTimeOffset'");
98  } else {
99  B2INFO("Exiting with failure");
100  return c_Failure;
101  }
102 
103 
104  /* 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency.
105  Same for all crystals. */
106  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::m_rf) * 1e3;
107 
108 
109 
110  // Set the rootfile name
111  bool minRunNumBool = false;
112  bool maxRunNumBool = false;
113  int minRunNum = -1;
114  int maxRunNum = -1;
115  int minExpNum = -1;
116  int maxExpNum = -1;
117  for (auto expRun : getRunList()) {
118  int expNumber = expRun.first;
119  int runNumber = expRun.second;
120  if (!minRunNumBool) {
121  minExpNum = expNumber;
122  minRunNum = runNumber;
123  minRunNumBool = true;
124  }
125  if (!maxRunNumBool) {
126  maxExpNum = expNumber;
127  maxRunNum = runNumber;
128  maxRunNumBool = true;
129  }
130  if (((minRunNum > runNumber) && (minExpNum >= expNumber)) ||
131  (minExpNum > expNumber)) {
132  minExpNum = expNumber;
133  minRunNum = runNumber;
134  }
135  if (((maxRunNum < runNumber) && (maxExpNum <= expNumber)) ||
136  (maxExpNum < expNumber))
137 
138  {
139  maxExpNum = expNumber;
140  maxRunNum = runNumber;
141  }
142  }
143 
144  string debugFilenameBase = "ECLCrystalTimeOffsetMerged";
145  B2INFO("debugFilenameBase = " << debugFilenameBase);
146  string runNumsString = string("_") + to_string(minExpNum) + "_" + to_string(minRunNum) + string("-") +
147  to_string(maxExpNum) + "_" + to_string(maxRunNum);
148  string debugFilename = debugFilenameBase + runNumsString + string(".root");
149 
150 
151 
152 
153  //------------------------------------------------------------------------
155  vector<float> bhabhaCalib = m_ECLCrystalTimeOffsetBhabha->getCalibVector();
156  vector<float> bhabhaCalibUnc = m_ECLCrystalTimeOffsetBhabha->getCalibUncVector(); // Negative uncertainties are flags of bad fits
157 
158  vector<float> bhabhaGammaCalib = m_ECLCrystalTimeOffsetBhabhaGamma->getCalibVector();
159  vector<float> bhabhaGammaCalibUnc = m_ECLCrystalTimeOffsetBhabhaGamma->getCalibUncVector(); // Negative uncertainties can be flags
160 
161  vector<float> cosmicCalib = m_ECLCrystalTimeOffsetCosmic->getCalibVector();
162  vector<float> cosmicCalibUnc = m_ECLCrystalTimeOffsetCosmic->getCalibUncVector();
163 
164  B2INFO("Loaded 'ECLCrystalTimeOffset*' calibrations");
165 
166  // Get the list of reference crystals as special case crystals
167  vector<short> refCrystalPerCrate = m_ECLReferenceCrystalPerCrateCalib->getReferenceCrystals();
168 
169  B2INFO("Loaded 'ECLReferenceCrystalPerCrateCalib' calibration");
170 
171  vector<float> crateCalib = m_ECLCrateTimeOffset->getCalibVector();
172 
173  B2INFO("Loaded 'ECLCrateTimeOffset' calibration");
174 
175 
176  // Get the previous crystal payloads
177  DBObjPtr<Belle2::ECLCrystalCalib> customPrevCrystalTimeObject("ECLCrystalTimeOffsetPreviousValues");
178  vector<float> prevValuesCrys(8736);
179  vector<float> prevValuesCrysUnc(8736);
180 
181  DBObjPtr<Belle2::ECLCrystalCalib> customPrevBhabhaCrystalTimeObject("ECLCrystalTimeOffsetBhabhaPreviousValues");
182  vector<float> prevBhabhaValuesCrys(8736);
183  vector<float> prevBhabhaValuesCrysUnc(8736);
184 
185  if (readPrevCrysPayload) {
186  //..Get vectors of values from the payloads
187  prevValuesCrys = customPrevCrystalTimeObject->getCalibVector();
188  prevValuesCrysUnc = customPrevCrystalTimeObject->getCalibUncVector();
189 
190  prevBhabhaValuesCrys = customPrevBhabhaCrystalTimeObject->getCalibVector();
191  prevBhabhaValuesCrysUnc = customPrevBhabhaCrystalTimeObject->getCalibUncVector();
192 
193  //..Print out a few values for quality control
194  B2INFO("Previous values read from database. Write out for their values for comparison");
195  for (int ic = 0; ic < 8736; ic += 500) {
196  B2INFO("ts custom previous payload: cellID " << ic + 1 << " " << prevValuesCrys[ic]);
197  }
198  B2INFO("Previous bhabha values read from database. Write out for their values for comparison");
199  for (int ic = 0; ic < 8736; ic += 500) {
200  B2INFO("ts custom previous bhabha payload: cellID " << ic + 1 << " " << prevBhabhaValuesCrys[ic]);
201  }
202  }
203 
204 
205 
206  //------------------------------------------------------------------------
234  //------------------------------------------------------------------------
235  /* Determine the bhabha and radiative bhabha calibration quality*/
236 
237  vector<bool> bhabhaCalibGoodQuality(m_numCrystals, false);
238  vector<bool> bhabhaGammaCalibGoodQuality(m_numCrystals, false);
239  for (int ic = 0; ic < m_numCrystals; ic++) {
240  /* Define a good bhabha calibration value. This is the
241  uncertainty on the calibration constant (e.g. fit mean),
242  not the guassian width of the timing distribution for
243  that crystal (i.e. resolution).
244 
245  Uncertainty stored as ticks so converted to ns for
246  more human readable cut.*/
247  if (bhabhaCalibUnc[ic] != 0
248  && fabs(bhabhaCalibUnc[ic])*TICKS_TO_NS < 2) {
249  // The limits on the uncertainties are not fully optimized !!!!!!!!!!!!!!!!!!!!!!!!!!
250  // Tighten the uncertainty cut after introduction of bhabha gamma calibration. Perhaps checking if the
251  // value is negative as that flags poor fits / the overall hist mean was used. In these cases,
252  // the bhabha gamma calibration should be used (if it has good stats) but until that calibration is
253  // properly made (not all zeros) then the bhabha calibration should be used instead as the only
254  // other option is the cosmic calibration.
255  bhabhaCalibGoodQuality[ic] = true ;
256  }
257  /* Define a good radiative bhabha calibration value, which
258  should be looser than the bhabha requirements since
259  the radiative bhabha calibration is the backup choice.
260  If the radiative bhabha calibration is not good then
261  we resort to using the cosmic calibration so the
262  radiative bhabha requirements should be fairly loose.
263 
264  Uncertainty stored as ticks so converted to ns for
265  more human readable cut.*/
266  if (bhabhaGammaCalibUnc[ic] != 0
267  && fabs(bhabhaGammaCalibUnc[ic])*TICKS_TO_NS <
268  3.5) { // The limits on the uncertainties are not fully optimized !!!!!!!!!!!!!!!!!!!!!!!!!!
269  bhabhaGammaCalibGoodQuality[ic] = true ;
270  }
271  }
272 
273 
274  /* Loop over all crystals and shift the cosmic crystal calibration by
275  the crate calibration because the cosmic calibration was made
276  assuming the crate calibration is 0.*/
277 
278  vector<float> crateCalibShort(m_numCrates);
279  vector<float> cosmicCalibShifted = cosmicCalib;
280  vector<float> cosmicCalibShiftedUnc = cosmicCalibUnc; // Do not modify uncertainty?
281  for (int ic = 0; ic < m_numCrystals; ic++) {
282  cosmicCalibShifted[ic] -= crateCalib[ic];
283  int crate_id_from_crystal = crystalMapper->getCrateID(ic + 1);
284  crateCalibShort[crate_id_from_crystal - 1] = crateCalib[ic];
285  }
286 
287  for (int icrate = 0; icrate < m_numCrates; icrate++) {
288  B2INFO("crate time calibration for crate " <<
289  icrate + 1 << " = " << crateCalibShort[icrate]);
290  }
291 
292 
293  /* Keep track of the mean value of the crystal calibration constants
294  for each crate. Exclude crystals with bad quality or that have
295  a value of 0 (except reference crystals). */
296  vector<float> meanGoodCrystalCalibPerCrate_bhabha(m_numCrates, 0);
297  vector<float> meanGoodCrystalCalibPerCrate_cosmic(m_numCrates, 0);
298 
299  /* Start num entries at 1 because we will skip all 0s but include
300  the reference crystal (ts=0) by just increasing the counter.
301  We only need a counter for the bhabha events since we want to
302  get the same sample of crystals for both the bhabha and cosmic
303  calibrations.*/
304  vector<short> numGoodCrystalCalibPerCrate_bhabha(m_numCrates, 1);
305 
306 
307  /* Loop over all crystals and save the crystal calibration information
308  for each crate. */
309  for (int ic_base1 = 1; ic_base1 <= m_numCrystals; ic_base1++) {
310  int crate_id_from_crystal = crystalMapper->getCrateID(ic_base1);
311 
312  if (!((bhabhaCalib[ic_base1 - 1] == 0) || (! bhabhaCalibGoodQuality[ic_base1 - 1]))) {
313  meanGoodCrystalCalibPerCrate_bhabha[crate_id_from_crystal - 1] += bhabhaCalib[ic_base1 - 1];
314  meanGoodCrystalCalibPerCrate_cosmic[crate_id_from_crystal - 1] += cosmicCalibShifted[ic_base1 - 1];
315  numGoodCrystalCalibPerCrate_bhabha[crate_id_from_crystal - 1]++ ;
316  }
317  }
318 
319  /* Complete the calculation of the mean by dividing the sum by the number of entries */
320  for (int icrate = 0; icrate < m_numCrates; icrate++) {
321  meanGoodCrystalCalibPerCrate_bhabha[icrate] /= numGoodCrystalCalibPerCrate_bhabha[icrate] ;
322  meanGoodCrystalCalibPerCrate_cosmic[icrate] /= numGoodCrystalCalibPerCrate_bhabha[icrate] ;
323  B2INFO("Mean bhabha crystal calibration for crate " << icrate + 1 << " = " << meanGoodCrystalCalibPerCrate_bhabha[icrate]);
324  B2INFO("Mean cosmic crystal calibration for crate " << icrate + 1 << " = " << meanGoodCrystalCalibPerCrate_cosmic[icrate]);
325  B2INFO("Number of good crystal calibrations for crate" << icrate + 1 << " = " << numGoodCrystalCalibPerCrate_bhabha[icrate]);
326  }
327 
328 
331  /* Loop over all crystals and correct the cosmic calibration by the cosmic calibration
332  for each crate by the average difference to the bhabha calibration. */
333  for (int ic_base1 = 1; ic_base1 <= m_numCrystals; ic_base1++) {
334  int crate_id_from_crystal = crystalMapper->getCrateID(ic_base1);
335  cosmicCalibShifted[ic_base1 - 1] += meanGoodCrystalCalibPerCrate_bhabha[crate_id_from_crystal - 1] -
336  meanGoodCrystalCalibPerCrate_cosmic[crate_id_from_crystal - 1] ;
337  }
338 
339 
340  //------------------------------------------------------------------------
350  vector<float> newCalib(m_numCrystals);
351  vector<float> newCalibUnc(m_numCrystals);
352 
353  vector<float> newBhabhaMinusCustomPrevCalib__cid(m_numCrystals);
354  vector<float> newBhabhaMinusCustomPrevBhabhaCalib__cid(m_numCrystals);
355  vector<float> newBhabhaMinusCustomPrevCalib;
356  vector<float> newBhabhaMinusCustomPrevBhabhaCalib;
357 
358 
359  // Loop over each crystal we want to check
360  for (int ic = 0; ic < m_numCrystals; ic++) {
361  B2DEBUG(29, "Merging crystal " << ic + 1);
362  int crate_id_from_crystal = crystalMapper->getCrateID(ic + 1);
363 
364  // Determine if the current crystal is a reference crystal
365  bool isRefCrys = false ;
366  for (int icrate = 0; icrate < m_numCrates; icrate++) {
367  if (ic + 1 == refCrystalPerCrate[icrate]) {
368  isRefCrys = true;
369  B2INFO("Crystal " << ic + 1 << " is a reference crystal");
370  break;
371  }
372  B2DEBUG(29, "Checked reference crystal for crate " << icrate + 1);
373  }
374  B2DEBUG(29, "Checked ref crystal");
375 
376  string whichCalibUsed = "bhabha";
377 
378  // By default set the new calibration to the bhabha calibration
379  newCalib[ic] = bhabhaCalib[ic];
380  newCalibUnc[ic] = fabs(bhabhaCalibUnc[ic]);
381 
382  double newMinusMerged ;
383  double newMinusBhabha ;
384 
385 
386  /* If the crystal is not a reference crystal and there is no bhabha
387  calib for that crystal then use a different calibration. If the
388  crystal is not a reference crystal and the crystal has a bad bhabha
389  crystal calibration uncertainty. */
390  if ((!isRefCrys) &&
391  ((bhabhaCalib[ic] == 0) || (! bhabhaCalibGoodQuality[ic]))) {
392 
393  /* If there is a radiative bhabha calibration for this crystal then
394  use it otherwise use the cosmic calibration */
395  if (bhabhaGammaCalib[ic] != 0 || bhabhaGammaCalibGoodQuality[ic]) {
396  newCalib[ic] = bhabhaGammaCalib[ic] ;
397  newCalibUnc[ic] = fabs(bhabhaGammaCalibUnc[ic]);
398  whichCalibUsed = "bhabhaGamma";
399  } else {
400  newCalib[ic] = cosmicCalibShifted[ic] ;
401  newCalibUnc[ic] = cosmicCalibShiftedUnc[ic];
402  whichCalibUsed = "cosmic";
403  }
404  } else {
405  newMinusMerged = (newCalib[ic] - prevValuesCrys[ic]) * TICKS_TO_NS;
406  newMinusBhabha = (newCalib[ic] - prevBhabhaValuesCrys[ic]) * TICKS_TO_NS;
407 
408  newBhabhaMinusCustomPrevCalib__cid[ic] = newMinusMerged;
409  newBhabhaMinusCustomPrevBhabhaCalib__cid[ic] = newMinusBhabha;
410  newBhabhaMinusCustomPrevCalib.push_back(newMinusMerged);
411  newBhabhaMinusCustomPrevBhabhaCalib.push_back(newMinusBhabha);
412  }
413 
414  newMinusMerged = (newCalib[ic] - prevValuesCrys[ic]) * TICKS_TO_NS;
415  newMinusBhabha = (newCalib[ic] - prevBhabhaValuesCrys[ic]) * TICKS_TO_NS;
416 
417 
418  B2INFO("Crystal " << ic + 1 << ", crate " << crate_id_from_crystal <<
419  "::: bhabha = " << bhabhaCalib[ic] << "+-" << bhabhaCalibUnc[ic] <<
420  " ticks, rad bhabha = " << bhabhaGammaCalib[ic] << "+-" << bhabhaGammaCalibUnc[ic] <<
421  " ticks, cosmic UNshifted = " << cosmicCalib[ic] << "+-" << cosmicCalibUnc[ic] <<
422  " ticks, cosmic shifted = " << cosmicCalibShifted[ic] << "+-" << cosmicCalibShiftedUnc[ic] <<
423  " ticks. ||| New calib = " << newCalib[ic] << " ticks, selected " << whichCalibUsed);
424  B2INFO(" Crystal " << ic + 1 << ", pre calib ts = " << prevValuesCrys[ic] <<
425  ", pre calib bhabha ts = " << prevBhabhaValuesCrys[ic] <<
426  " ticks ||| new - pre calib ts = " << newMinusMerged <<
427  ", new - pre calib bhabha ts = " << newMinusBhabha << " ns");
428  } // end of setting new calibration for each crystal
429 
430 
431 
432  //------------------------------------------------------------------------
433  //..Histograms of calibrations
434 
435  // Just in case, we remember the current TDirectory so we can return to it
436  TDirectory* executeDir = gDirectory;
437 
438  TString fname = debugFilename;
439  TFile hfile(fname, "recreate");
440 
441  TString htitle = "ECLCrystalTimeOffsetBhabha;cellID;Bhabha ts values [ticks]";
442  TH1F* bhabhaPayload = new TH1F("bhabhaPayload", htitle, m_numCrystals, 1, 8737);
443 
444  htitle = "ECLCrystalTimeOffsetBhabhaGamma;cellID;Radiative bhabha ts values [ticks]";
445  TH1F* bhabhaGammaPayload = new TH1F("bhabhaGammaPayload", htitle, m_numCrystals, 1, 8737);
446 
447  htitle = "ECLCrystalTimeOffsetCosmic : unshifted;cellID;Unshifted cosmic ts values [ticks]";
448  TH1F* cosmicPayload = new TH1F("cosmicUnshiftedPayload", htitle, m_numCrystals, 1, 8737);
449 
450  htitle = "ECLCrystalTimeOffsetCosmic : shifted;cellID;Shifted cosmic ts values [ticks]";
451  TH1F* cosmicShiftedPayload = new TH1F("cosmicShiftedPayload", htitle, m_numCrystals, 1, 8737);
452 
453  htitle = "New ECLCrystalTimeOffset;cellID;New (merged) ts values [ticks]";
454  TH1F* newPayload = new TH1F("newPayload", htitle, m_numCrystals, 1, 8737);
455 
456 
457  htitle = "ECLCrystalTimeOffsetBhabha : 'pre-calib' values;cellID;Pre-calibration bhabha ts values [ticks]";
458  TH1F* tsCustomPrevBhabha_payload = new TH1F("tsCustomPrevBhabha_payload", htitle, m_numCrystals, 1, 8737);
459 
460  htitle = "ECLCrystalTimeOffset : 'pre-calib' values;cellID;Pre-calibration merged ts values [ticks]";
461  TH1F* tsCustomPrev_payload = new TH1F("tsCustomPrev_payload", htitle, m_numCrystals, 1, 8737);
462 
463 
464 
465  htitle = "-Only for crystals where the new bhabha ts value is used-;cellID;ts(new | bhabha) - ts(old = 'pre-calib' | merged) [ns]";
466  TH1F* newBhabhaMinusCustomPrev__cid = new TH1F("newBhabhaMinusCustomPrev__cid", htitle, m_numCrystals, 1, 8737);
467 
468  htitle = "-Only for crystals where the new bhabha ts value is used-;cellID;ts(new | bhabha) - ts(old = 'pre-calib' | bhabha) [ns]";
469  TH1F* newBhabhaMinusCustomPrevBhabha__cid = new TH1F("newBhabhaMinusCustomPrevBhabha__cid", htitle, m_numCrystals, 1, 8737);
470 
471  htitle = "-Only for crystals where the new bhabha ts value is used-;ts(new | bhabha) - ts(old = 'pre-calib' | merged) [ns];Number of crystals";
472  auto tsNewBhabha_MINUS_tsCustomPrev = new TH1F("TsNewBhabha_MINUS_TsCustomPrev", htitle, 285, -69.5801, 69.5801);
473 
474  htitle = "-Only for crystals where the new bhabha ts value is used-;ts(new | bhabha) - ts(old = 'pre-calib' | bhabha) [ns];Number of crystals";
475  auto tsNewBhabha_MINUS_tsCustomPrevBhabha = new TH1F("TsNewBhabha_MINUS_TsCustomPrevBhabha", htitle, 285, -69.5801, 69.5801);
476 
477 
478  for (int cellID = 1; cellID <= m_numCrystals; cellID++) {
479  bhabhaPayload->SetBinContent(cellID, bhabhaCalib[cellID - 1]);
480  bhabhaPayload->SetBinError(cellID, fabs(bhabhaCalibUnc[cellID - 1]));
481 
482  bhabhaGammaPayload->SetBinContent(cellID, bhabhaGammaCalib[cellID - 1]);
483  bhabhaGammaPayload->SetBinError(cellID, fabs(bhabhaGammaCalibUnc[cellID - 1]));
484 
485  cosmicPayload->SetBinContent(cellID, cosmicCalib[cellID - 1]);
486  cosmicPayload->SetBinError(cellID, cosmicCalibUnc[cellID - 1]);
487 
488  cosmicShiftedPayload->SetBinContent(cellID, cosmicCalibShifted[cellID - 1]);
489  cosmicShiftedPayload->SetBinError(cellID, cosmicCalibShiftedUnc[cellID - 1]);
490 
491  newPayload->SetBinContent(cellID, newCalib[cellID - 1]);
492  newPayload->SetBinError(cellID, newCalibUnc[cellID - 1]);
493 
494  // Comparisons of change in ts
495  if (readPrevCrysPayload) {
496  newBhabhaMinusCustomPrev__cid->SetBinContent(cellID, newBhabhaMinusCustomPrevCalib__cid[cellID - 1]);
497  newBhabhaMinusCustomPrev__cid->SetBinError(cellID, 0);
498 
499  newBhabhaMinusCustomPrevBhabha__cid->SetBinContent(cellID, newBhabhaMinusCustomPrevBhabhaCalib__cid[cellID - 1]);
500  newBhabhaMinusCustomPrevBhabha__cid->SetBinError(cellID, 0);
501 
502 
503  tsCustomPrev_payload->SetBinContent(cellID, prevValuesCrys[cellID - 1]);
504  tsCustomPrev_payload->SetBinError(cellID, prevValuesCrysUnc[cellID - 1]);
505 
506  tsCustomPrevBhabha_payload->SetBinContent(cellID, prevBhabhaValuesCrys[cellID - 1]);
507  tsCustomPrevBhabha_payload->SetBinError(cellID, fabs(prevBhabhaValuesCrysUnc[cellID - 1]));
508  }
509  }
510 
511  if (readPrevCrysPayload) {
512  for (long unsigned int ib = 0; ib < newBhabhaMinusCustomPrevCalib.size(); ib++) {
513  tsNewBhabha_MINUS_tsCustomPrev->Fill(newBhabhaMinusCustomPrevCalib[ib]);
514  }
515 
516  for (long unsigned int ib = 0; ib < newBhabhaMinusCustomPrevBhabhaCalib.size(); ib++) {
517  tsNewBhabha_MINUS_tsCustomPrevBhabha->Fill(newBhabhaMinusCustomPrevBhabhaCalib[ib]);
518  }
519  }
520 
521 
522  hfile.cd();
523  hfile.Write();
524  hfile.Close();
525  B2INFO("Debugging histograms written to " << fname);
526  // Go back to original TDirectory
527  executeDir->cd();
528 
529 
530  /* Re-save the new bhabha calibrations to a payload. The bhabha calibrations have NOT
531  been changed; however, the best way of controling which payloads get get uploaded to
532  the GT as part of the prommp calibrations is to save a second copy in this calibration
533  directory since we don't want to save all the payloads from the previous directory.
534  See the "<calibration>.save_payloads = <True/False>" code in the airflow steering file.*/
535  ECLCrystalCalib* bhabhaCalibCopy = new ECLCrystalCalib();
536  bhabhaCalibCopy->setCalibVector(bhabhaCalib, bhabhaCalibUnc);
537  saveCalibration(bhabhaCalibCopy, "ECLCrystalTimeOffsetBhabha");
538 
539  // Save the new merged calibrations to a payload
540  ECLCrystalCalib* newCrystalTimes = new ECLCrystalCalib();
541  newCrystalTimes->setCalibVector(newCalib, newCalibUnc);
542  saveCalibration(newCrystalTimes, "ECLCrystalTimeOffset");
543  return c_OK;
544 }
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.
@ c_Failure
Failed =3 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
General DB object to store one calibration number per ECL crystal.
const std::vector< float > & getCalibUncVector() const
Get vector of uncertainties on calibration constants.
const std::vector< float > & getCalibVector() const
Get vector of calibration constants.
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
static constexpr double m_rf
accelerating RF, http://ptep.oxfordjournals.org/content/2013/3/03A006.full.pdf
DBObjPtr< ECLCrystalCalib > m_ECLCrystalTimeOffsetBhabha
ECLCrystalTimeOffsetBhabha payload that we want to read from the DB.
DBObjPtr< ECLReferenceCrystalPerCrateCalib > m_ECLReferenceCrystalPerCrateCalib
ECLReferenceCrystalPerCrateCalib payload that we want to read from the DB.
static constexpr int m_numCrates
Number of Crates expected.
DBObjPtr< ECLCrystalCalib > m_ECLCrateTimeOffset
ECLCrateTimeOffset payload that we want to read from the DB.
static constexpr int m_numCrystals
Number of Crystals expected.
bool readPrevCrysPayload
< Read the previous crystal payload values for comparison
DBObjPtr< ECLCrystalCalib > m_ECLCrystalTimeOffsetCosmic
ECLCrystalTimeOffsetCosmic payload that we want to read from the DB.
DBObjPtr< ECLCrystalCalib > m_ECLCrystalTimeOffsetBhabhaGamma
ECLCrystalTimeOffsetBhabhaGamma payload that we want to read from the DB.
Abstract base class for different kinds of events.
Struct containing exp number and run number.
Definition: Splitter.h:51