Belle II Software  release-06-02-00
eclBhabhaTAlgorithm.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 #include <ecl/calibration/eclBhabhaTAlgorithm.h>
10 #include <ecl/dbobjects/ECLCrystalCalib.h>
11 #include <ecl/dbobjects/ECLReferenceCrystalPerCrateCalib.h>
12 #include <ecl/digitization/EclConfiguration.h>
13 #include <ecl/utility/ECLChannelMapper.h>
14 
15 #include <framework/database/DBObjPtr.h>
16 #include <framework/database/DBStore.h>
17 #include <framework/datastore/StoreObjPtr.h>
18 #include <framework/datastore/DataStore.h>
19 #include <framework/dataobjects/EventMetaData.h>
20 
21 #include "TH2F.h"
22 #include "TFile.h"
23 #include "TF1.h"
24 #include "TROOT.h"
25 
26 using namespace std;
27 using namespace Belle2;
28 using namespace ECL;
29 
30 eclBhabhaTAlgorithm::eclBhabhaTAlgorithm():
31  // Parameters
32  CalibrationAlgorithm("ECLBhabhaTCollector"),
33  cellIDLo(1),
34  cellIDHi(8736),
35  meanCleanRebinFactor(1),
36  meanCleanCutMinFactor(0),
37  crateIDLo(1),
38  crateIDHi(52),
39  savePrevCrysPayload(false),
40  readPrevCrysPayload(false),
41  debugOutput(true),
42  debugFilenameBase("eclBhabhaTAlgorithm"), // base of filename (without ".root")
43  collectorName("ECLBhabhaTCollector"),
44  refCrysPerCrate{ -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
45  -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
46  -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
47  -1, -1, -1, -1, -1, -1, -1} //,
48 {
49  setDescription("Calculate time offsets from bhabha events by fitting gaussian function to the (t - T0) difference.");
50 }
51 
52 
53 
55 {
57  gROOT->SetBatch();
58 
59 
61  B2INFO("eclBhabhaTAlgorithm parameters:");
62  B2INFO("cellIDLo = " << cellIDLo);
63  B2INFO("cellIDHi = " << cellIDHi);
64  B2INFO("meanCleanRebinFactor = " << meanCleanRebinFactor);
65  B2INFO("meanCleanCutMinFactor = " << meanCleanCutMinFactor);
66  B2INFO("crateIDLo = " << crateIDLo);
67  B2INFO("crateIDHi = " << crateIDHi);
68  B2INFO("savePrevCrysPayload = " << savePrevCrysPayload);
69  B2INFO("readPrevCrysPayload = " << readPrevCrysPayload);
70  B2INFO("refCrysPerCrate = {");
71  for (int crateTest = 0; crateTest < 52 - 1; crateTest++) {
72  B2INFO(refCrysPerCrate[crateTest] << ",");
73  }
74  B2INFO(refCrysPerCrate[52 - 1] << "}");
75 
76 
77  /* Histogram with the data collected by eclBhabhaTCollectorModule */
78 
79  auto TimevsCrysPrevCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibPrevCrystCalib");
80  auto TimevsCratePrevCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCratePrevCrateCalibPrevCrystCalib");
81  auto TimevsCrysNoCalibrations = getObjectPtr<TH2F>("TimevsCrysNoCalibrations");
82  auto TimevsCrysPrevCrateCalibNoCrystCalib = getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibNoCrystCalib");
83  auto TimevsCrateNoCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCrateNoCrateCalibPrevCrystCalib");
84 
85  // Collect other plots just for reference - combines all the runs for these plots.
86  auto cutflow = getObjectPtr<TH1F>("cutflow");
87 
88 
89  // Define new plots to make
90  // New calibration constant values minus older values from previous iteration plotted as a function of the crystal or crate id
91  unique_ptr<TH1F> tsNew_MINUS_tsOld__cid(new TH1F("TsNew_MINUS_TsOld__cid",
92  ";cell id; ts(new|bhabha) - ts(previous iteration|merged) [ns]", 8736,
93  1, 8736 + 1));
94  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld__crateID(new TH1F("tcrateNew_MINUS_tcrateOld__crateID",
95  ";crate id; tcrate(new | bhabha) - tcrate(previous iteration | merged) [ns]",
96  52, 1, 52 + 1));
97  unique_ptr<TH1F> tsNew_MINUS_tsCustomPrev__cid(new TH1F("TsNew_MINUS_TsCustomPrev__cid",
98  ";cell id; ts(new|bhabha) - ts(old = 'before 1st iter'|merged) [ns]",
99  8736, 1, 8736 + 1));
100  unique_ptr<TH1F> tsNew_MINUS_tsOldBhabha__cid(new TH1F("TsNew_MINUS_TsOldBhabha__cid",
101  ";cell id; ts(new|bhabha) - ts(previous iteration|bhabha) [ns]", 8736,
102  1, 8736 + 1));
103 
104 
105  // Histogram of the new time constant values minus values from previous iteration
106  unique_ptr<TH1F> tsNew_MINUS_tsOld(new TH1F("TsNew_MINUS_TsOld",
107  ";ts(new | bhabha) - ts(previous iteration | merged) [ns];Number of crystals",
108  201, -10.05, 10.05));
109  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld(new TH1F("tcrateNew_MINUS_tcrateOld",
110  ";tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
111  201, -10.05, 10.05));
112  unique_ptr<TH1F> tsNew_MINUS_tsCustomPrev(new TH1F("TsNew_MINUS_TsCustomPrev",
113  ";ts(new | bhabha) - ts(old = 'before 1st iter' | merged) [ns];Number of crystals",
114  285, -69.5801, 69.5801));
115  unique_ptr<TH1F> tsNew_MINUS_tsOldBhabha(new TH1F("TsNew_MINUS_TsOldBhabha",
116  ";ts(new | bhabha) - ts(previous iteration | bhabha) [ns];Number of crystals",
117  201, -10.05, 10.05));
118 
119 
120 
121  // Histograms just for crates and collecting for all the different runs rather than run by run.
122  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld_allRuns(new TH1F("tcrateNew_MINUS_tcrateOld_allRuns",
123  "Crate time constant changes over all runs : tcrate(new) uncertainty < 0.1ns;tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
124  201, -10.05, 10.05));
125 
126  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld_allRuns_allCrates(new TH1F("tcrateNew_MINUS_tcrateOld_allRuns_allCrates",
127  "Crate time constant changes over all runs : all crates;tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
128  201, -10.05, 10.05));
129 
130  unique_ptr<TH1I> num_tcrates_perRun(new TH1I("num_tcrates_perRun",
131  "Number of good tcrates in each run;Run number;Number of good tcrates",
132  6000, 0, 6000));
133 
134  unique_ptr<TH2F> tcrateNew_MINUS_tcrateOld__vs__runNum(new TH2F("tcrateNew_MINUS_tcrateOld__vs__runNum",
135  "Crate time constant changes vs run number : tcrate(new) uncertainty < 0.1ns;Run number;tcrate(new) - tcrate(previous iteration) [ns]",
136  6000, 0, 6000, 21, -10.5, 10.5));
137 
138 
139 
140 
141 
142  if (!TimevsCrysNoCalibrations) return c_Failure;
143 
146  TFile* histfile = 0;
147  TFile* histExtraCrateInfofile = 0; // Only created if needed for crates
148 
149 
150  unique_ptr<TTree> tree_crystal(new TTree("tree_crystal", "Debug data from bhabha time calibration algorithm for crystals"));
151 
152  unique_ptr<TTree> tree_crate(new TTree("tree_crate", "Debug data from bhabha time calibration algorithm for crates"));
153  int tree_cid;
154 
155  // Vector of time offsets to be saved in the database.
156  vector<float> t_offsets;
157  // Vector of time offset uncertainties to be saved in the database.
158  vector<float> t_offsets_unc;
159  vector<float> t_offsets_prev; // previous time offsets
160 
161 
162  int minNumEntries = 40;
163  int minNumEntriesCrateConvergence = 1000;
164 
165 
166  double mean = 0;
167  double sigma = -1;
168  double mean_unc = 0;
169  int crystalCalibSaved = 0;
170  double tsPrev = 0;
171 
172 
173  bool minRunNumBool = false;
174  bool maxRunNumBool = false;
175  int minRunNum = -1;
176  int maxRunNum = -1;
177  int minExpNum = -1;
178  int maxExpNum = -1;
179  for (auto expRun : getRunList()) {
180  int expNumber = expRun.first;
181  int runNumber = expRun.second;
182  if (!minRunNumBool) {
183  minExpNum = expNumber;
184  minRunNum = runNumber;
185  minRunNumBool = true;
186  }
187  if (!maxRunNumBool) {
188  maxExpNum = expNumber;
189  maxRunNum = runNumber;
190  maxRunNumBool = true;
191  }
192  if (((minRunNum > runNumber) && (minExpNum >= expNumber)) ||
193  (minExpNum > expNumber)) {
194  minExpNum = expNumber;
195  minRunNum = runNumber;
196  }
197  if (((maxRunNum < runNumber) && (maxExpNum <= expNumber)) ||
198  (maxExpNum < expNumber))
199 
200  {
201  maxExpNum = expNumber;
202  maxRunNum = runNumber;
203  }
204  }
205 
206  B2INFO("debugFilenameBase = " << debugFilenameBase);
207  string runNumsString = string("_") + to_string(minExpNum) + "_" + to_string(minRunNum) + string("-") +
208  to_string(maxExpNum) + "_" + to_string(maxRunNum);
209  string debugFilename = debugFilenameBase + runNumsString + string(".root");
210  string extraCratedebugFilename = debugFilenameBase + string("_cratesAllRuns.root");
211 
212 
213  // Need to load information about the event/run/experiment to get the right database information
214  // Will be used for:
215  // * ECLChannelMapper (to map crystal to crates)
216  // * crystal payload updating for iterating crystal and crate fits
217  int eventNumberForCrates = 1;
218 
220  // simulate the initialize() phase where we can register objects in the DataStore
222  evtPtr.registerInDataStore();
224  // now construct the event metadata
225  evtPtr.construct(eventNumberForCrates, minRunNum, minExpNum);
226  // and update the database contents
227  DBStore& dbstore = DBStore::Instance();
228  dbstore.update();
229  // this is only needed it the payload might be intra-run dependent,
230  // that is if it might change during one run as well
231  dbstore.updateEvent();
232 
233 
234  B2INFO("Uploading payload for exp " << minExpNum << ", run " << minRunNum << ", event " << eventNumberForCrates);
235  updateDBObjPtrs(eventNumberForCrates, minRunNum, minExpNum);
236  unique_ptr<ECLChannelMapper> crystalMapper(new ECL::ECLChannelMapper());
237  crystalMapper->initFromDB();
238 
239 
240  //------------------------------------------------------------------------
241  //..Read payloads from database
242  B2INFO("Reading payloads: ECLCrystalTimeOffset and ECLCrateTimeOffset");
243  DBObjPtr<Belle2::ECLCrystalCalib> crystalTimeObject("ECLCrystalTimeOffset");
244  DBObjPtr<Belle2::ECLCrystalCalib> crateTimeObject("ECLCrateTimeOffset");
245 
246  //..Get vectors of values from the payloads
247  vector<float> currentValuesCrys = crystalTimeObject->getCalibVector();
248  vector<float> currentUncCrys = crystalTimeObject->getCalibUncVector();
249  vector<float> currentValuesCrate = crateTimeObject->getCalibVector();
250  vector<float> currentUncCrate = crateTimeObject->getCalibUncVector();
251 
252  //..Print out a few values for quality control
253  B2INFO("Values read from database. Write out for their values for comparison against those from tcol");
254  for (int ic = 0; ic < 8736; ic += 500) {
255  B2INFO("ts: cellID " << ic + 1 << " " << currentValuesCrys[ic] << " +/- " << currentUncCrys[ic]);
256  B2INFO("tcrate: cellID " << ic + 1 << " " << currentValuesCrate[ic] << " +/- " << currentUncCrate[ic]);
257  }
258 
259 
260  //..Read in the previous crystal payload values
261  vector<float> prevValuesCrys(8736);
262  if (readPrevCrysPayload) {
263  DBObjPtr<Belle2::ECLCrystalCalib> customPrevCrystalTimeObject("ECLCrystalTimeOffsetPreviousValues");
264  //..Get vectors of values from the payloads
265  prevValuesCrys = customPrevCrystalTimeObject->getCalibVector();
266 
267  //..Print out a few values for quality control
268  B2INFO("Previous values read from database. Write out for their values for comparison against those from tcol");
269  for (int ic = 0; ic < 8736; ic += 500) {
270  B2INFO("ts custom previous payload: cellID " << ic + 1 << " " << prevValuesCrys[ic]);
271  }
272  }
273 
274 
275  //..Read bhabha payloads from database
276  B2INFO("Reading payloads: ECLCrystalTimeOffsetBhabha");
277  DBObjPtr<Belle2::ECLCrystalCalib> crystalBhabhaTimeObject("ECLCrystalTimeOffsetBhabha");
278 
279  //..Get vectors of values from the payloads
280  vector<float> currentBhabhaValuesCrys = crystalBhabhaTimeObject->getCalibVector();
281  vector<float> currentBhabhaUncCrys = crystalBhabhaTimeObject->getCalibUncVector();
282 
283  //..Print out a few values for quality control
284  for (int ic = 0; ic < 8736; ic += 500) {
285  B2INFO("ts bhabha: cellID " << ic + 1 << " " << currentBhabhaValuesCrys[ic] << " +/- " << currentBhabhaUncCrys[ic]);
286  }
287 
288 
289  //------------------------------------------------------------------------
290  //..Get the reference crystal information
291  auto refCrysIDzeroingCrate = getObjectPtr<TH1F>("refCrysIDzeroingCrate");
292 
293  // crystal index for the crystals (one per crate) that is used as the reference crystal. This one has
294  // ts defined as zero. The crystal id runs 1...8636, not starting at 0.
295  B2INFO("Extract reference crystals from collector histogram.");
296  vector <short> crystalIDreferenceUntested;
297  for (int bin = 1; bin <= 8736; bin++) {
298  if (refCrysIDzeroingCrate->GetBinContent(bin) > 0.5) {
299  crystalIDreferenceUntested.push_back(bin);
300  }
301  }
302 
303  // Output for the user the crystal id to be used as a reference
304  // and also state which crate the crystal is in.
305  B2INFO("Reference crystals to define as having ts=0. Base 1 counting for both crates and crystals");
306  for (long unsigned int crysRefCounter = 0; crysRefCounter < crystalIDreferenceUntested.size(); crysRefCounter++) {
307  int crys_id = crystalIDreferenceUntested[crysRefCounter] ;
308  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
309  B2INFO(" crystal " << crys_id << " is a reference for crate " << crate_id_from_crystal);
310  }
311 
312  // Check that the reference crystals make sense. There should be exactly one per crate but
313  // since the crystal calibration executes over multiple runs, it is possible that the
314  // reference crystals could change but we don't want to allow this.
315  B2INFO("Checking number of reference crystals");
316  B2INFO("Number of reference crystals = " << crystalIDreferenceUntested.size());
317 
318  // Check that there are exactly 52 reference crystals
319  if (crystalIDreferenceUntested.size() != 52) {
320  B2FATAL("The number of reference crystals does not equal 52, which is one per crate");
321  return c_Failure;
322  } else {
323  B2INFO("Number of reference crystals is 52 as required");
324  }
325 
326  // Count the number of reference crystals for each crate to make sure that there is exactly
327  // one reference crystal for each crate as defined by the payload/database.
328  // also fill the final vector that maps the crate id to the reference crystal id
329  vector <short> crateIDsNumRefCrystalsUntested(52, 0);
330  vector <short> crystalIDReferenceForZeroTs(52, 0);
331 
332  for (long unsigned int crysRefCounter = 0; crysRefCounter < crystalIDreferenceUntested.size(); crysRefCounter++) {
333  int crys_id = crystalIDreferenceUntested[crysRefCounter] ;
334  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
335  crateIDsNumRefCrystalsUntested[crate_id_from_crystal - 1]++;
336  crystalIDReferenceForZeroTs[crate_id_from_crystal - 1] = crys_id;
337  }
338  B2INFO("crystalIDReferenceForZeroTs = {");
339  for (int crateTest = 0; crateTest < 52 - 1; crateTest++) {
340  B2INFO(crystalIDReferenceForZeroTs[crateTest] << ",");
341  }
342  B2INFO(crystalIDReferenceForZeroTs[52 - 1] << "}");
343 
344 
345  // Make sure that there is only one reference crystal per crate as defined by the payload/database
346  for (int crateTest = 0; crateTest < 52; crateTest++) {
347  if (crateIDsNumRefCrystalsUntested[crateTest] != 1) {
348  B2FATAL("Crate " << crateTest + 1 << " (base 1) has " << crateIDsNumRefCrystalsUntested[crateTest] << " reference crystals");
349  return c_Failure;
350  }
351  }
352  B2INFO("All reference crystals are reasonably mapped one crystal to one crate for all crates");
353 
354 
355 
356  B2INFO("Extract reference crystals from algorithm steering script if provided. If user inputs custom values via steering script for this algorithm, they are only applied after all the tests are performed on the values from the histogram and override the histogram valuees. User can adjust just a single crystal if desired. Use -1 to indicate that a crystal is not to be modified. Position of crystal in list determines the crate to which the crystal is meant to be associated.");
357 
358  /* Test if the user wants to modify the reference crystals. This will probably be done very rarely,
359  perhaps less than once per year. If the user does want to change one or more reference
360  crystals then perform the checks again to make sure that there is still just one reference crystal
361  per crate after the modifications to the payload values with the user values.*/
362  bool userSetRefCrysPerCrate = false ;
363  for (int crateTest = 0; crateTest < 52; crateTest++) {
364  if (refCrysPerCrate[crateTest] != -1) {
365  crystalIDReferenceForZeroTs[crateTest] = refCrysPerCrate[crateTest] ;
366  B2INFO("Crate " << crateTest + 1 << " (base 1) new reference crystal = " << crystalIDReferenceForZeroTs[crateTest]);
367  userSetRefCrysPerCrate = true ;
368  }
369  }
370  if (userSetRefCrysPerCrate) {
371  B2INFO("User changed reference crystals via steering script");
372 
373  // Validate crystals per crate again with the new user set values
374  fill(crateIDsNumRefCrystalsUntested.begin(), crateIDsNumRefCrystalsUntested.end(), 0);
375  for (long unsigned int crysRefCounter = 0; crysRefCounter < 52; crysRefCounter++) {
376  int crys_id = crystalIDReferenceForZeroTs[crysRefCounter] ;
377  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
378  crateIDsNumRefCrystalsUntested[crate_id_from_crystal - 1]++;
379  }
380  for (int crateTest = 0; crateTest < 52; crateTest++) {
381  if (crateIDsNumRefCrystalsUntested[crateTest] != 1) {
382  B2FATAL("Crate " << crateTest + 1 << " (base 1) has " << crateIDsNumRefCrystalsUntested[crateTest] << " reference crystals");
383  return c_Failure;
384  }
385  }
386  B2INFO("All reference crystals are reasonably mapped one crystal to one crate for all crates after changes made by user steering script.");
387 
388  // Save the information to the payload if there is at least one crate
389  // reference crystal that has been modified by the user steering file
391  refCrystalsCalib->setCalibVector(crystalIDReferenceForZeroTs);
392  saveCalibration(refCrystalsCalib, "ECLReferenceCrystalPerCrateCalib");
393  B2INFO("Created reference crystal per crate payload: ECLReferenceCrystalPerCrateCalib");
394  } else {
395  B2INFO("User did not change reference crystals via steering script");
396  }
397 
398 
399  //------------------------------------------------------------------------
400  //..Start looking at timing information
401 
402  B2INFO("Debug output rootfile: " << debugFilename);
403  histfile = new TFile(debugFilename.c_str(), "recreate");
404 
405 
406  TimevsCrysPrevCrateCalibPrevCrystCalib ->Write();
407  TimevsCratePrevCrateCalibPrevCrystCalib->Write();
408  TimevsCrysNoCalibrations ->Write();
409  TimevsCrysPrevCrateCalibNoCrystCalib ->Write();
410  TimevsCrateNoCrateCalibPrevCrystCalib ->Write();
411 
412  cutflow->Write();
413 
414 
415  if (debugOutput) {
416  tree_crystal->Branch("cid", &tree_cid)->SetTitle("Cell ID, 1..8736");
417  tree_crystal->Branch("ts", &mean)->SetTitle("Time offset mean, ts, ns");
418  tree_crystal->Branch("tsUnc", &mean_unc)->SetTitle("Error of time ts mean, ns.");
419  tree_crystal->Branch("tsSigma", &sigma)->SetTitle("Sigma of time ts distribution, ns");
420  tree_crystal->Branch("crystalCalibSaved",
421  &crystalCalibSaved)->SetTitle("0=crystal skipped, 1=crystal calib saved (num entries based)");
422  tree_crystal->Branch("tsPrev", &tsPrev)->SetTitle("Previous crystal time offset, ts, ns");
423  tree_crystal->SetAutoSave(10);
424  }
425 
426 
427  double hist_tmin = TimevsCrysNoCalibrations->GetYaxis()->GetXmin();
428  double hist_tmax = TimevsCrysNoCalibrations->GetYaxis()->GetXmax();
429 
430  double time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
431  double time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
432 
433  B2INFO("hist_tmin = " << hist_tmin);
434  B2INFO("hist_tmax = " << hist_tmax);
435 
436  /* 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency.
437  Same for all crystals. */
438  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::m_rf) * 1e3;
439 
440  // The ts and tcrate database values are filled once per tcol instance so count the number of times that the database values
441  // were summed together by the histogram merging process and extract out the original values again.
442  auto databaseCounter = getObjectPtr<TH1I>("databaseCounter");
443  float numTimesFilled = databaseCounter->GetBinContent(1);
444  B2INFO("Number of times database histograms were merged = " << numTimesFilled);
445 
446 
447  auto TsDatabase = getObjectPtr<TH1F>("TsDatabase");
448  auto TsDatabaseUnc = getObjectPtr<TH1F>("TsDatabaseUnc");
449  for (int i = 1; i <= 8736; i++) {
450  t_offsets.push_back(TsDatabase->GetBinContent(i) / numTimesFilled);
451  t_offsets_prev.push_back(TsDatabase->GetBinContent(i) / numTimesFilled);
452 
453  B2INFO("t_offsets_prev (last iter) at crysID " << i << " = " << t_offsets_prev[i - 1]);
454 
455  t_offsets_unc.push_back(TsDatabaseUnc->GetBinContent(i) / numTimesFilled);
456  }
457 
458 
459  /* CRYSTAL CORRECTIONS */
460 
461  /* Make a 1D histogram of the number of hits per crystal. This will help with the validations
462  to make sure that all the crystals had enough hits and to look for problems.*/
463  TH1D* h_crysHits = TimevsCrysPrevCrateCalibNoCrystCalib->ProjectionX("h_crysHits");
464  h_crysHits->SetTitle("Hits per crystal;Crystal id");
465 
466  histfile->WriteTObject(h_crysHits, "h_crysHits");
467 
468 
469  // Loop over all the crystals for doing the crystal calibation
470  for (int crys_id = cellIDLo; crys_id <= cellIDHi; crys_id++) {
471  crystalCalibSaved = 0;
472 
473  double database_mean = 0;
474  double database_mean_unc = 0;
475 
476  B2INFO("Crystal id = " << crys_id);
477 
478 
479  vector<bool> ts_new_was_set(8736, false);
480 
481 
482  /* Determining which bins to mask out for mean calculation
483  */
484 
485  TH1D* h_time = TimevsCrysPrevCrateCalibNoCrystCalib->ProjectionY((string("h_time_psi__") + to_string(crys_id)).c_str(),
486  crys_id, crys_id);
487  TH1D* h_timeMask = (TH1D*)h_time->Clone();
488  TH1D* h_timeMasked = (TH1D*)h_time->Clone((string("h_time_psi_masked__") + to_string(crys_id)).c_str());
489  TH1D* h_timeRebin = (TH1D*)h_time->Clone();
490 
491  // Do rebinning and cleaning of some bins but only if the user selection values call for it since it slows the code down
492  if (meanCleanRebinFactor != 1 || meanCleanCutMinFactor != 1) {
493 
494  h_timeRebin->Rebin(meanCleanRebinFactor);
495 
496  h_timeMask->Scale(0.0); // set all bins to being masked off
497 
498  time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
499  time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
500 
501  // Find value of bin with max value
502  double histRebin_max = h_timeRebin->GetMaximum();
503 
504  bool maskedOutNonZeroBin = false;
505  // Loop over all bins to find those with content less than a certain threshold. Mask the non-rebinned histogram for the corresponding bins
506  for (int bin = 1; bin <= h_timeRebin->GetNbinsX(); bin++) {
507  for (int rebinCounter = 1; rebinCounter <= meanCleanRebinFactor; rebinCounter++) {
508  int nonRebinnedBinNumber = (bin - 1) * meanCleanRebinFactor + rebinCounter;
509  if (nonRebinnedBinNumber < h_time->GetNbinsX()) {
510  if (h_timeRebin->GetBinContent(bin) >= histRebin_max * meanCleanCutMinFactor) {
511  h_timeMask->SetBinContent(nonRebinnedBinNumber, 1);
512 
513  // Save the lower and upper edges of the rebin histogram time range for fitting purposes
514  double x_lower = h_timeRebin->GetXaxis()->GetBinLowEdge(bin);
515  double x_upper = h_timeRebin->GetXaxis()->GetBinUpEdge(bin);
516  if (x_lower < time_fit_min) {
517  time_fit_min = x_lower;
518  }
519  if (x_upper > time_fit_max) {
520  time_fit_max = x_upper;
521  }
522 
523  } else {
524  if (h_time->GetBinContent(nonRebinnedBinNumber) > 0) {
525  B2DEBUG(22, "Setting bin " << nonRebinnedBinNumber << " from " << h_timeMasked->GetBinContent(nonRebinnedBinNumber) << " to 0");
526  maskedOutNonZeroBin = true;
527  }
528  h_timeMasked->SetBinContent(nonRebinnedBinNumber, 0);
529  }
530  }
531  }
532  }
533  B2INFO("Bins with non-zero values have been masked out: " << maskedOutNonZeroBin);
534  h_timeMasked->ResetStats();
535  h_timeMask->ResetStats();
536 
537  }
538 
539  // Calculate mean from masked histogram
540  double default_meanMasked = h_timeMasked->GetMean();
541  //double default_meanMasked_unc = h_timeMasked->GetMeanError();
542  B2INFO("default_meanMasked = " << default_meanMasked);
543 
544 
545  // Get the overall mean and standard deviation of the distribution within the plot. This doesn't require a fit.
546  double default_mean = h_time->GetMean();
547  double default_mean_unc = h_time->GetMeanError();
548  double default_sigma = h_time->GetStdDev();
549 
550  B2INFO("Fitting crystal between " << time_fit_min << " and " << time_fit_max);
551 
552  // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
553  TF1* gaus = new TF1("func", "gaus(0)", time_fit_min, time_fit_max);
554  gaus->SetParNames("numCrystalHitsNormalization", "mean", "sigma");
555  /*
556  gaus->ReleaseParameter(0); // number of crystals
557  gaus->ReleaseParameter(1); // mean
558  gaus->ReleaseParameter(2); // standard deviation
559  */
560 
561  double hist_max = h_time->GetMaximum();
562 
563  //=== Estimate initial value of sigma as std dev.
564  double stddev = h_time->GetStdDev();
565  sigma = stddev;
566  mean = default_mean;
567 
568  //=== Setting parameters for initial iteration
569  gaus->SetParameter(0, hist_max / 2.);
570  gaus->SetParameter(1, mean);
571  gaus->SetParameter(2, sigma);
572  // L -- Use log likelihood method
573  // I -- Use integral of function in bin instead of value at bin center // not using
574  // R -- Use the range specified in the function range
575  // B -- Fix one or more parameters with predefined function // not using
576  // Q -- Quiet mode
577 
578  h_timeMasked->Fit(gaus, "LQR"); // L for likelihood, R for x-range, Q for fit quiet mode
579 
580  double fit_mean = gaus->GetParameter(1);
581  double fit_mean_unc = gaus->GetParError(1);
582  double fit_sigma = gaus->GetParameter(2);
583 
584  double meanDiff = fit_mean - default_mean;
585  double meanUncDiff = fit_mean_unc - default_mean_unc;
586  double sigmaDiff = fit_sigma - default_sigma;
587 
588  bool good_fit = false;
589 
590  if ((fabs(meanDiff) > 10) ||
591  (fabs(meanUncDiff) > 10) ||
592  (fabs(sigmaDiff) > 10) ||
593  (fit_mean_unc > 0.09) ||
594  (fit_sigma < 0.1) ||
595  (fit_mean < time_fit_min) ||
596  (fit_mean > time_fit_max)) {
597  B2INFO("Crystal id = " << crys_id);
598  B2INFO("fit mean, default mean = " << fit_mean << ", " << default_mean);
599  B2INFO("fit mean unc, default mean unc = " << fit_mean_unc << ", " << default_mean_unc);
600  B2INFO("fit sigma, default sigma = " << fit_sigma << ", " << default_sigma);
601 
602  B2INFO("crystal fit mean - hist mean = " << meanDiff);
603  B2INFO("fit mean unc. - hist mean unc. = " << meanUncDiff);
604  B2INFO("fit sigma - hist sigma = " << sigmaDiff);
605 
606  B2INFO("fit_mean = " << fit_mean);
607  B2INFO("time_fit_min = " << time_fit_min);
608  B2INFO("time_fit_max = " << time_fit_max);
609 
610  if (fabs(meanDiff) > 10) B2INFO("fit mean diff too large");
611  if (fabs(meanUncDiff) > 10) B2INFO("fit mean unc diff too large");
612  if (fabs(sigmaDiff) > 10) B2INFO("fit mean sigma diff too large");
613  if (fit_mean_unc > 0.09) B2INFO("fit mean unc too large");
614  if (fit_sigma < 0.1) B2INFO("fit sigma too small");
615 
616  } else {
617  good_fit = true;
618  }
619 
620 
621 
622  // Set the tree_crystal values - ignore fit values !!!!!!!!!!!!!!!!
623  sigma = default_sigma;
624 
625 
626  int numEntries = h_time->GetEntries();
627  /* If number of entries in histogram is greater than X then use the statistical information from the data otherwise
628  leave crystal uncalibrated. Histograms are still shown though. ALSO require the that fits are good.*/
629  if ((numEntries >= minNumEntries) && good_fit) {
630  crystalCalibSaved = 1;
631  database_mean = fit_mean;
632  database_mean_unc = fit_mean_unc;
633  } else {
634  database_mean = default_mean;
635  database_mean_unc = -fabs(default_mean_unc);
636  }
637 
638  if (numEntries < minNumEntries) B2INFO("Number of entries less than minimum");
639  if (numEntries == 0) B2INFO("Number of entries == 0");
640 
641 
642  tree_cid = crys_id;
643 
644  // For the database, convert back from ns to ADC ticks.
645  t_offsets[crys_id - 1] = database_mean / TICKS_TO_NS;
646  t_offsets_unc[crys_id - 1] = database_mean_unc / TICKS_TO_NS;
647 
648 
649  histfile->WriteTObject(h_time, (string("h_time_psi") + to_string(crys_id)).c_str());
650  histfile->WriteTObject(h_timeMasked, (string("h_time_psi_masked") + to_string(crys_id)).c_str());
651 
652  mean = database_mean;
653  mean_unc = database_mean_unc;
654 
655  tsPrev = t_offsets_prev[crys_id - 1] * TICKS_TO_NS;
656 
657  delete gaus;
658  tree_crystal->Fill();
659  }
660 
661 
662  // Shift the crystal time calibration constants by the reference crystal calibration constant values
663  if (cellIDLo <= cellIDHi) {
664  vector <double> tsRefCID ;
665  B2INFO("crystal times before shift");
666  for (int crate_id = 1; crate_id <= 52; crate_id++) {
667  tsRefCID.push_back(t_offsets[ crystalIDReferenceForZeroTs[crate_id - 1] - 1 ]);
668  B2INFO("crystal time [crystal = " << crystalIDReferenceForZeroTs[crate_id - 1] << ", crate = " << crate_id << " (base 1)] = " <<
669  t_offsets[ crystalIDReferenceForZeroTs[crate_id - 1] - 1 ] << " ticks");
670  }
671 
672  B2INFO("crystal times after shift wrt reference crystal");
673  for (int crys_id = 1; crys_id <= 8736; crys_id++) {
674  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
675  B2INFO("crystal time before shift [crystal = " << crys_id << ", crate = " << crate_id_from_crystal << " (base 1)] = " <<
676  t_offsets[crys_id - 1] << " +- " << t_offsets_unc[crys_id - 1] << " ticks");
677 
678  /* Shift the crystal time constant by that of the reference crystal, but only if
679  there were values to shift. If there were no entries, ts=0 and ts_unc=0, which
680  are special values so do not shift these crystals. */
681  if (t_offsets[crys_id - 1] == 0 && t_offsets_unc[crys_id - 1] == 0) {
682  B2INFO("crystal time after shift [crystal = " << crys_id << ", crate = " << crate_id_from_crystal << " (base 1)] = " <<
683  t_offsets[crys_id - 1] << " ticks. No change because ts=0 and ts_unc=0 (no entries).");
684  } else {
685  t_offsets[crys_id - 1] = t_offsets[crys_id - 1] - tsRefCID[crate_id_from_crystal - 1];
686  B2INFO("crystal time after shift [crystal = " << crys_id << ", crate = " << crate_id_from_crystal << " (base 1)] = " <<
687  t_offsets[crys_id - 1] << " ticks");
688  }
689 
690 
691  // Fill histograms with the difference in the ts values between iterations
692  double tsDiff_ns = (t_offsets[crys_id - 1] - t_offsets_prev[crys_id - 1]) * TICKS_TO_NS;
693  double tsDiffBhabha_ns = -999;
694  if (readPrevCrysPayload) {
695  tsDiffBhabha_ns = (t_offsets[crys_id - 1] - currentBhabhaValuesCrys[crys_id - 1]) * TICKS_TO_NS;
696  }
697 
698  B2INFO("Crystal " << crys_id << ": ts new bhabha - old merged = (" <<
699  t_offsets[crys_id - 1] << " - " << t_offsets_prev[crys_id - 1] <<
700  ") ticks * " << TICKS_TO_NS << " ns/tick = " << tsDiff_ns << " ns");
701  B2INFO("Crystal " << crys_id << ": ts new bhabha - old bhabha = (" <<
702  t_offsets[crys_id - 1] << " - " << currentBhabhaValuesCrys[crys_id - 1] <<
703  ") ticks * " << TICKS_TO_NS << " ns/tick = " << tsDiffBhabha_ns << " ns");
704 
705  tsNew_MINUS_tsOld__cid->SetBinContent(crys_id, tsDiff_ns);
706  tsNew_MINUS_tsOld__cid->SetBinError(crys_id, 0);
707  tsNew_MINUS_tsOld__cid->ResetStats();
708 
709  tsNew_MINUS_tsOld->Fill(tsDiff_ns);
710 
711 
712  tsNew_MINUS_tsOldBhabha__cid->SetBinContent(crys_id, tsDiffBhabha_ns);
713  tsNew_MINUS_tsOldBhabha__cid->SetBinError(crys_id, 0);
714  tsNew_MINUS_tsOldBhabha__cid->ResetStats();
715 
716  tsNew_MINUS_tsOldBhabha->Fill(tsDiffBhabha_ns);
717 
718 
719 
720  /* Fill histograms with the difference in the ts values from this iteration
721  and the previous values read in from the payload. */
722  double tsDiffCustomOld_ns = -999;
723  if (readPrevCrysPayload) {
724  tsDiffCustomOld_ns = (t_offsets[crys_id - 1] - prevValuesCrys[crys_id - 1]) * TICKS_TO_NS;
725  B2INFO("Crystal " << crys_id << ": ts new bhabha - 'before 1st iter' merged = (" <<
726  t_offsets[crys_id - 1] << " - " << prevValuesCrys[crys_id - 1] <<
727  ") ticks * " << TICKS_TO_NS << " ns/tick = " << tsDiffCustomOld_ns << " ns");
728  }
729  tsNew_MINUS_tsCustomPrev__cid->SetBinContent(crys_id, tsDiffCustomOld_ns);
730  tsNew_MINUS_tsCustomPrev__cid->SetBinError(crys_id, 0);
731  tsNew_MINUS_tsCustomPrev__cid->ResetStats();
732 
733  tsNew_MINUS_tsCustomPrev->Fill(tsDiffCustomOld_ns);
734 
735  }
736 
737  // Save the histograms to the output root file
738  histfile->WriteTObject(tsNew_MINUS_tsOld__cid.get(), "tsNew_MINUS_tsOld__cid");
739  histfile->WriteTObject(tsNew_MINUS_tsOld.get(), "tsNew_MINUS_tsOld");
740 
741  histfile->WriteTObject(tsNew_MINUS_tsCustomPrev__cid.get(), "tsNew_MINUS_tsCustomPrev__cid");
742  histfile->WriteTObject(tsNew_MINUS_tsCustomPrev.get(), "tsNew_MINUS_tsCustomPrev");
743 
744  histfile->WriteTObject(tsNew_MINUS_tsOldBhabha__cid.get(), "tsNew_MINUS_tsOldBhabha__cid");
745  histfile->WriteTObject(tsNew_MINUS_tsOldBhabha.get(), "tsNew_MINUS_tsOldBhabha");
746  }
747 
748 
749  //..Store previous crystal calibration constants to payload under different
750  // names so that they can be read in for comparison later. These are temporary
751  // payloads that are only used for plotting purposes while running the calibration
752  // and does not need to be added to any global tag.
753  if (savePrevCrysPayload) {
754  ECLCrystalCalib* crysTCalib_prev = new ECLCrystalCalib();
755  crysTCalib_prev->setCalibVector(currentValuesCrys, currentUncCrys);
756 
757  ECLCrystalCalib* crysBhabhaTCalib_prev = new ECLCrystalCalib();
758  crysBhabhaTCalib_prev->setCalibVector(currentBhabhaValuesCrys, currentBhabhaUncCrys);
759 
760 
761  // Save the information to the payload if there is at least one crystal
762  // begin calibrated.
763  if (cellIDLo <= cellIDHi) {
764  saveCalibration(crysTCalib_prev, "ECLCrystalTimeOffsetPreviousValues");
765  B2INFO("Previous overall crystal payload made");
766 
767  saveCalibration(crysBhabhaTCalib_prev, "ECLCrystalTimeOffsetBhabhaPreviousValues");
768  B2INFO("Previous bhabha crystal payload made");
769  }
770  }
771 
772 
773  ECLCrystalCalib* BhabhaTCalib = new ECLCrystalCalib();
774  BhabhaTCalib->setCalibVector(t_offsets, t_offsets_unc);
775 
776  // Save the information to the payload if there is at least one crystal
777  // begin calibrated.
778  if (cellIDLo <= cellIDHi) {
779  saveCalibration(BhabhaTCalib, "ECLCrystalTimeOffset");
780  saveCalibration(BhabhaTCalib, "ECLCrystalTimeOffsetBhabha");
781  B2DEBUG(22, "crystal payload made");
782  }
783 
784 
785  B2DEBUG(22, "end of crystal start of crate corrections .....");
786 
787 
788  //==============================================================
789  /* CRATE CORRECTIONS */
790 
791  hist_tmin = TimevsCrateNoCrateCalibPrevCrystCalib->GetYaxis()->GetXmin();
792  hist_tmax = TimevsCrateNoCrateCalibPrevCrystCalib->GetYaxis()->GetXmax();
793 
794  B2DEBUG(22, "Found min/max of X axis of TimevsCrateNoCrateCalibPrevCrystCalib");
795 
796  // Vector of time offsets to be saved in the database.
797 
798  auto TcrateDatabase = getObjectPtr<TH1F>("TcrateDatabase");
799 
800 
801  B2DEBUG(22, "Retrieved Ts and Tcrate histograms from tcol root file");
802 
803 
804  vector<float> tcrate_mean_new(52, 0.0);
805  vector<float> tcrate_mean_unc_new(52, 0.0);
806  vector<float> tcrate_sigma_new(52, 0.0);
807  vector<float> tcrate_mean_prev(52, 0.0);
808  // vector<float> tcrate_sigma_prev(52, 0.0); // currently not used
809  vector<bool> tcrate_new_was_set(52, false);
810  vector<bool> tcrate_new_goodQuality(52, false);
811 
812  B2DEBUG(22, "crate vectors set");
813 
814 
815 
816  // Crate time calibration constants are saved per crystal so read them per crystal
817  // and save as one entry per crate in the array
818  for (int crys_id = 1; crys_id <= 8736; crys_id++) {
819  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
820  tcrate_mean_prev[crate_id_from_crystal - 1] = TcrateDatabase->GetBinContent(crys_id) / numTimesFilled;
821  }
822 
823 
824  B2INFO("Print out previous crate time calibration constants to make sure they match from the two different sources.");
825  for (int crate_id = 1; crate_id <= 52; crate_id++) {
826  B2INFO("tcrate_mean_prev[crate " << crate_id << " (base 1)] = " << tcrate_mean_prev[crate_id - 1]);
827 
828  int thisRefCellID = crystalIDReferenceForZeroTs[crate_id - 1];
829  B2INFO("tcrate from payload: ref cellID " << thisRefCellID << " " << currentValuesCrate[thisRefCellID - 1] << " +/- " <<
830  currentUncCrate[thisRefCellID - 1]);
831  }
832 
833 
834  /* Read in the histogram about the crate calibration constant differences between iterations
835  if it exists. This way the histograms can be updated after each run.*/
836  if (crateIDLo <= crateIDHi) {
837  TFile* histExtraCrateInfofile_dummy = 0; // Only created if needed for crates
838  B2INFO("Debug output rootfile used for crate iterations: " << extraCratedebugFilename);
839  histExtraCrateInfofile_dummy = new TFile(extraCratedebugFilename.c_str(), "UPDATE");
840 
841 
842  // If the histogram already exists then read it in and recreate the file so that we can save an updated version of the histogram.
843  TKey* key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld_allRuns");
844  if (key != 0) {
845  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld_allRuns");
846  tcrateNew_MINUS_tcrateOld_allRuns->Add(h);
847  }
848 
849  key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
850  if (key != 0) {
851  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
852  tcrateNew_MINUS_tcrateOld_allRuns_allCrates->Add(h);
853  }
854 
855  key = histExtraCrateInfofile_dummy->FindKey("num_tcrates_perRun");
856  if (key != 0) {
857  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("num_tcrates_perRun");
858  num_tcrates_perRun->Add(h);
859  }
860 
861  key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld__vs__runNum");
862  if (key != 0) {
863  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld__vs__runNum");
864  tcrateNew_MINUS_tcrateOld__vs__runNum->Add(h);
865  }
866 
867  histExtraCrateInfofile_dummy->Close();
868 
869  /* After reading in all the histograms, recreate the root file from empty so that the
870  histograms can be made again with the updated values.*/
871  histExtraCrateInfofile = new TFile(extraCratedebugFilename.c_str(), "recreate");
872  }
873 
874 
875 
876  for (int crate_id = crateIDLo; crate_id <= crateIDHi; crate_id++) {
877 
878  B2DEBUG(22, "Start of crate id = " << crate_id);
879 
880  TH1D* h_time_crate = TimevsCrateNoCrateCalibPrevCrystCalib->ProjectionY("h_time_psi_crate", crate_id, crate_id);
881  TH1D* h_time_crate_mask = (TH1D*)h_time_crate->Clone();
882  TH1D* h_time_crate_masked = (TH1D*)h_time_crate->Clone();
883  TH1D* h_time_crate_rebin = (TH1D*)h_time_crate->Clone();
884 
885 
886  // Do rebinning and cleaning of some bins but only if the user selection values call for it since it slows the code down
887  if (meanCleanRebinFactor != 1 || meanCleanCutMinFactor != 1) {
888 
889  h_time_crate_rebin->Rebin(meanCleanRebinFactor);
890  h_time_crate_mask->Scale(0.0); // set all bins to being masked off
891 
892  time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
893  time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
894 
895  // Find value of bin with max value
896  double histRebin_max = h_time_crate_rebin->GetMaximum();
897 
898  bool maskedOutNonZeroBin = false;
899  // Loop over all bins to find those with content less than a certain threshold. Mask the non-rebinned histogram for the corresponding bins
900  for (int bin = 1; bin <= h_time_crate_rebin->GetNbinsX(); bin++) {
901  for (int rebinCounter = 1; rebinCounter <= meanCleanRebinFactor; rebinCounter++) {
902  int nonRebinnedBinNumber = (bin - 1) * meanCleanRebinFactor + rebinCounter;
903  if (nonRebinnedBinNumber < h_time_crate->GetNbinsX()) {
904  if (h_time_crate_rebin->GetBinContent(bin) >= histRebin_max * meanCleanCutMinFactor) {
905  h_time_crate_mask->SetBinContent(nonRebinnedBinNumber, 1);
906 
907  // Save the lower and upper edges of the rebin histogram time range for fitting purposes
908  double x_lower = h_time_crate_rebin->GetXaxis()->GetBinLowEdge(bin);
909  double x_upper = h_time_crate_rebin->GetXaxis()->GetBinUpEdge(bin);
910  if (x_lower < time_fit_min) {
911  time_fit_min = x_lower;
912  }
913  if (x_upper > time_fit_max) {
914  time_fit_max = x_upper;
915  }
916  } else {
917  if (h_time_crate->GetBinContent(nonRebinnedBinNumber) > 0) {
918  B2DEBUG(22, "Setting bin " << nonRebinnedBinNumber << " from " << h_time_crate_masked->GetBinContent(
919  nonRebinnedBinNumber) << " to 0");
920  maskedOutNonZeroBin = true;
921  }
922  h_time_crate_masked->SetBinContent(nonRebinnedBinNumber, 0);
923  }
924  }
925  }
926  }
927  B2INFO("Bins with non-zero values have been masked out: " << maskedOutNonZeroBin);
928  h_time_crate_masked->ResetStats();
929  h_time_crate_mask->ResetStats();
930 
931  }
932 
933 
934 
935  B2DEBUG(22, "crate loop - projected h_time_psi_crate");
936 
937 
938  double default_mean_crate = h_time_crate_masked->GetMean();
939  double default_mean_crate_unc = h_time_crate_masked->GetMeanError();
940  double default_sigma_crate = h_time_crate_masked->GetStdDev();
941  B2INFO("Fitting crate between " << time_fit_min << " and " << time_fit_max);
942  TF1* gaus = new TF1("func", "gaus(0)", time_fit_min, time_fit_max);
943  gaus->SetParNames("numCrateHisNormalization", "mean", "sigma");
944  double hist_max = h_time_crate->GetMaximum();
945  double stddev = h_time_crate->GetStdDev();
946  double sigma_crate = stddev;
947  double mean_crate = default_mean_crate;
948  gaus->SetParameter(0, hist_max / 2.);
949  gaus->SetParameter(1, mean_crate);
950  gaus->SetParameter(2, sigma_crate);
951 
952  h_time_crate_masked->Fit(gaus, "LQR"); // L for likelihood, R for x-range, Q for fit quiet mode
953 
954  double fit_mean_crate = gaus->GetParameter(1);
955  double fit_mean_crate_unc = gaus->GetParError(1);
956  double fit_sigma_crate = gaus->GetParameter(2);
957 
958  double meanDiff = fit_mean_crate - default_mean_crate;
959  double meanUncDiff = fit_mean_crate_unc - default_mean_crate_unc;
960  double sigmaDiff = fit_sigma_crate - default_sigma_crate;
961 
962  bool good_fit = false;
963 
964  B2DEBUG(22, "Crate id = " << crate_id << " with crate mean = " << default_mean_crate << " +- " << fit_mean_crate_unc);
965 
966  if ((fabs(meanDiff) > 7) ||
967  (fabs(meanUncDiff) > 7) ||
968  (fabs(sigmaDiff) > 7) ||
969  (fit_mean_crate_unc > 3) ||
970  (fit_sigma_crate < 0.1) ||
971  (fit_mean_crate < time_fit_min) ||
972  (fit_mean_crate > time_fit_max)) {
973  B2INFO("Crate id = " << crate_id);
974  B2INFO("fit mean, default mean = " << fit_mean_crate << ", " << default_mean_crate);
975  B2INFO("fit sigma, default sigma = " << fit_sigma_crate << ", " << default_sigma_crate);
976 
977  B2INFO("crate fit mean - hist mean = " << meanDiff);
978  B2INFO("fit mean unc. - hist mean unc. = " << meanUncDiff);
979  B2INFO("fit sigma - hist sigma = " << sigmaDiff);
980  B2INFO("fit_mean_crate = " << fit_mean_crate);
981  B2INFO("time_fit_min = " << time_fit_min);
982  B2INFO("time_fit_max = " << time_fit_max);
983  } else {
984  good_fit = true;
985  }
986 
987  int numEntries = h_time_crate->GetEntries();
988  B2INFO("Number entries = " << numEntries);
989  double database_mean_crate = 0;
990  double database_mean_crate_unc = 0;
991  if ((numEntries >= minNumEntries) && good_fit) {
992  database_mean_crate = fit_mean_crate;
993  database_mean_crate_unc = fit_mean_crate_unc;
994 
995  if ((numEntries >= minNumEntriesCrateConvergence) && (fit_mean_crate_unc < 0.1)) {
996  tcrate_new_goodQuality[crate_id - 1] = true;
997  }
998  } else {
999  database_mean_crate = default_mean_crate;
1000  database_mean_crate_unc = default_mean_crate_unc;
1001  }
1002  if (numEntries == 0) {
1003  database_mean_crate = 0;
1004  database_mean_crate_unc = 0;
1005  }
1006 
1007  tcrate_mean_new[crate_id - 1] = database_mean_crate;
1008  tcrate_mean_unc_new[crate_id - 1] = database_mean_crate_unc;
1009  tcrate_sigma_new[crate_id - 1] = fit_sigma_crate;
1010  tcrate_new_was_set[crate_id - 1] = true;
1011 
1012 
1013  histfile->WriteTObject(h_time_crate, (string("h_time_psi_crate") + to_string(crate_id)).c_str());
1014  histfile->WriteTObject(h_time_crate_masked, (string("h_time_psi_crate_masked") + to_string(crate_id)).c_str());
1015  histfile->WriteTObject(h_time_crate_rebin, (string("h_time_psi_crate_rebinned") + to_string(crate_id)).c_str());
1016 
1017  delete gaus;
1018  }
1019 
1020  B2DEBUG(22, "crate histograms made");
1021 
1022 
1023  // Save database for crates
1024  // Vector of time offsets to be saved in the database.
1025  vector<float> t_offsets_crate;
1026  // Vector of time offset uncertainties to be saved in the database.
1027  vector<float> t_offsets_crate_unc;
1028  for (int i = 1; i <= 8736; i++) {
1029  t_offsets_crate.push_back(0);
1030  t_offsets_crate_unc.push_back(0);
1031  }
1032 
1033 
1034  for (int crys_id = 1; crys_id <= 8736; crys_id++) {
1035  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
1036  if (tcrate_new_was_set[crate_id_from_crystal - 1]) {
1037  t_offsets_crate[crys_id - 1] = tcrate_mean_new[crate_id_from_crystal - 1] / TICKS_TO_NS;
1038  t_offsets_crate_unc[crys_id - 1] = tcrate_mean_unc_new[crate_id_from_crystal - 1] / TICKS_TO_NS;
1039 
1040  } else {
1041  t_offsets_crate[crys_id - 1] = tcrate_mean_prev[crate_id_from_crystal - 1];
1042  B2INFO("used old crate mean but zeroed uncertainty since not saved in root file");
1043  }
1044  }
1045 
1046 
1047  // Fill histograms with the difference in the tcrate values
1048  if (crateIDLo <= crateIDHi) {
1049  for (int crate_id = crateIDLo; crate_id <= crateIDHi; crate_id++) {
1050  // tcrate_mean_new already in ns, but tcrate_mean_prev in ticks.
1051  double tCrateDiff_ns = tcrate_mean_new[crate_id - 1] - (tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS);
1052  B2INFO("Crate " << crate_id << ": tcrate new - previous iteration = "
1053  << tcrate_mean_new[crate_id - 1]
1054  << " - " << tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS
1055  << " = " << tCrateDiff_ns << " ns");
1056  tcrateNew_MINUS_tcrateOld__crateID->SetBinContent(crate_id, tCrateDiff_ns);
1057  tcrateNew_MINUS_tcrateOld__crateID->SetBinError(crate_id, 0);
1058  tcrateNew_MINUS_tcrateOld__crateID->ResetStats();
1059 
1060  tcrateNew_MINUS_tcrateOld->Fill(tCrateDiff_ns);
1061 
1062  // Save the histograms monitoring the change between iterations
1063  tcrateNew_MINUS_tcrateOld_allRuns_allCrates->Fill(tCrateDiff_ns);
1064  if (tcrate_new_goodQuality[crate_id - 1]) {
1065  tcrateNew_MINUS_tcrateOld_allRuns->Fill(tCrateDiff_ns);
1066  num_tcrates_perRun->Fill(minRunNum);
1067  tcrateNew_MINUS_tcrateOld__vs__runNum->Fill(minRunNum, tCrateDiff_ns);
1068  }
1069  }
1070 
1071  // Save the histograms to the output root file
1072  histfile->WriteTObject(tcrateNew_MINUS_tcrateOld__crateID.get(), "tcrateNew_MINUS_tcrateOld__crateID");
1073  histfile->WriteTObject(tcrateNew_MINUS_tcrateOld.get(), "tcrateNew_MINUS_tcrateOld");
1074 
1075  // Save the histograms ot the crate iterations file
1076  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld_allRuns.get(), "tcrateNew_MINUS_tcrateOld_allRuns");
1077  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld_allRuns_allCrates.get(),
1078  "tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
1079  histExtraCrateInfofile->WriteTObject(num_tcrates_perRun.get(), "num_tcrates_perRun");
1080  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld__vs__runNum.get(), "tcrateNew_MINUS_tcrateOld__vs__runNum");
1081  }
1082 
1083 
1084  ECLCrystalCalib* BhabhaTCrateCalib = new ECLCrystalCalib();
1085  BhabhaTCrateCalib->setCalibVector(t_offsets_crate, t_offsets_crate_unc);
1086 
1087 
1088  // Save the information to the payload if there is at least one crate
1089  // begin calibrated.
1090  if (crateIDLo <= crateIDHi) {
1091  saveCalibration(BhabhaTCrateCalib, "ECLCrateTimeOffset");
1092  B2DEBUG(22, "crate payload made");
1093 
1094  histExtraCrateInfofile->Close();
1095  }
1096 
1097 
1098  int tree_crateid;
1099  int tree_runNum;
1100  double tree_tcrate_mean;
1101  double tree_tcrate_mean_unc;
1102  double tree_tcrate_sigma;
1103  double tree_tcrate_meanPrev;
1104 
1105  tree_crate->Branch("runNum", &tree_runNum)->SetTitle("Run number, 0..infinity and beyond!");
1106  tree_crate->Branch("crateid", &tree_crateid)->SetTitle("Crate id, 1..52");
1107  tree_crate->Branch("tcrate", &tree_tcrate_mean)->SetTitle("Crate time offset mean, tcrate, ns");
1108  tree_crate->Branch("tcratePrev", &tree_tcrate_meanPrev)->SetTitle("Previous crate time offset mean, tcrate, ns");
1109  tree_crate->Branch("tcrate_unc", &tree_tcrate_mean_unc)->SetTitle("Error of time tcrate mean, ns.");
1110  tree_crate->Branch("tcrate_sigma", &tree_tcrate_sigma)->SetTitle("Sigma of time tcrate distribution, ns");
1111  tree_crate->SetAutoSave(10);
1112 
1113 
1114  for (auto expRun : getRunList()) {
1115  // Key command to make sure your DBObjPtrs are correct
1116  B2INFO("run num, exp num: " << expRun.second << ", " << expRun.first);
1117  int runNumber = expRun.second;
1118 
1119  for (int crate_id = 1; crate_id <= 52; crate_id++) {
1120  if (tcrate_new_was_set[crate_id - 1]) {
1121  tree_runNum = runNumber;
1122  tree_crateid = crate_id;
1123  tree_tcrate_mean = tcrate_mean_new[crate_id - 1];
1124  tree_tcrate_mean_unc = tcrate_mean_unc_new[crate_id - 1];
1125  tree_tcrate_sigma = tcrate_sigma_new[crate_id - 1];
1126  tree_tcrate_meanPrev = tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS;
1127  tree_crate->Fill();
1128  }
1129  }
1130  }
1131 
1132  B2DEBUG(22, "end of crate corrections .....");
1133 
1134  tree_crystal->Write();
1135  tree_crate->Write();
1136 
1137  histfile->Close();
1138 
1139  B2INFO("Finished talgorithm");
1140  return c_OK;
1141 }
1142 
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.
Singleton class to cache database objects.
Definition: DBStore.h:32
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:52
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:92
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.
General DB object to store one reference crystal per per ECL crate for calibration purposes.
void setCalibVector(const std::vector< short > &refCrystals)
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
int crateIDHi
Fit crates with crateID0 in the inclusive range [crateIDLo,crateIDHi].
int cellIDHi
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
int cellIDLo
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
bool debugOutput
Save every histogram and fitted function to debugFilename.
double meanCleanRebinFactor
Rebinning factor for mean calculation.
int crateIDLo
Fit crates with crateID0 in the inclusive range [crateIDLo,crateIDHi].
double meanCleanCutMinFactor
After rebinning, create a mask for bins that have values less than meanCleanCutMinFactor times the ma...
bool savePrevCrysPayload
Save the previous crystal payload values for comparison.
bool readPrevCrysPayload
Read the previous crystal payload values for comparison.
EResult calibrate() override
..Run algorithm on events
int refCrysPerCrate[52]
List of crystals, one per crate, used as reference time for crystal time calibration.
std::string debugFilenameBase
Name of file with debug output, eclBhabhaTAlgorithm.root by default.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:118
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:26
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:140
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:77
Abstract base class for different kinds of events.