Belle II Software  release-05-02-19
eclBhabhaTAlgorithm.cc
1 #include <ecl/calibration/eclBhabhaTAlgorithm.h>
2 #include <ecl/dbobjects/ECLCrystalCalib.h>
3 #include <ecl/dbobjects/ECLReferenceCrystalPerCrateCalib.h>
4 #include <ecl/digitization/EclConfiguration.h>
5 #include <ecl/utility/ECLChannelMapper.h>
6 
7 #include <framework/database/DBObjPtr.h>
8 #include <framework/database/DBStore.h>
9 #include <framework/datastore/StoreObjPtr.h>
10 #include <framework/datastore/DataStore.h>
11 #include <framework/dataobjects/EventMetaData.h>
12 
13 #include "TH2F.h"
14 #include "TFile.h"
15 #include "TF1.h"
16 #include "TROOT.h"
17 
18 using namespace std;
19 using namespace Belle2;
20 using namespace ECL;
21 
22 eclBhabhaTAlgorithm::eclBhabhaTAlgorithm():
23  // Parameters
24  CalibrationAlgorithm("ECLBhabhaTCollector"),
25  cellIDLo(1),
26  cellIDHi(8736),
27  meanCleanRebinFactor(1),
28  meanCleanCutMinFactor(0),
29  crateIDLo(1),
30  crateIDHi(52),
31  savePrevCrysPayload(false),
32  readPrevCrysPayload(false),
33  debugOutput(true),
34  debugFilenameBase("eclBhabhaTAlgorithm"), // base of filename (without ".root")
35  collectorName("ECLBhabhaTCollector"),
36  refCrysPerCrate{ -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
37  -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
38  -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
39  -1, -1, -1, -1, -1, -1, -1} //,
40 {
41  setDescription("Calculate time offsets from bhabha events by fitting gaussian function to the (t - T0) difference.");
42 }
43 
44 
45 
47 {
49  gROOT->SetBatch();
50 
51 
53  B2INFO("eclBhabhaTAlgorithm parameters:");
54  B2INFO("cellIDLo = " << cellIDLo);
55  B2INFO("cellIDHi = " << cellIDHi);
56  B2INFO("meanCleanRebinFactor = " << meanCleanRebinFactor);
57  B2INFO("meanCleanCutMinFactor = " << meanCleanCutMinFactor);
58  B2INFO("crateIDLo = " << crateIDLo);
59  B2INFO("crateIDHi = " << crateIDHi);
60  B2INFO("savePrevCrysPayload = " << savePrevCrysPayload);
61  B2INFO("readPrevCrysPayload = " << readPrevCrysPayload);
62  B2INFO("refCrysPerCrate = {");
63  for (int crateTest = 0; crateTest < 52 - 1; crateTest++) {
64  B2INFO(refCrysPerCrate[crateTest] << ",");
65  }
66  B2INFO(refCrysPerCrate[52 - 1] << "}");
67 
68 
69  /* Histogram with the data collected by eclBhabhaTCollectorModule */
70 
71  auto TimevsCrysPrevCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibPrevCrystCalib");
72  auto TimevsCratePrevCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCratePrevCrateCalibPrevCrystCalib");
73  auto TimevsCrysNoCalibrations = getObjectPtr<TH2F>("TimevsCrysNoCalibrations");
74  auto TimevsCrysPrevCrateCalibNoCrystCalib = getObjectPtr<TH2F>("TimevsCrysPrevCrateCalibNoCrystCalib");
75  auto TimevsCrateNoCrateCalibPrevCrystCalib = getObjectPtr<TH2F>("TimevsCrateNoCrateCalibPrevCrystCalib");
76 
77  // Collect other plots just for reference - combines all the runs for these plots.
78  auto cutflow = getObjectPtr<TH1F>("cutflow");
79 
80 
81  // Define new plots to make
82  // New calibration constant values minus older values from previous iteration plotted as a function of the crystal or crate id
83  unique_ptr<TH1F> tsNew_MINUS_tsOld__cid(new TH1F("TsNew_MINUS_TsOld__cid",
84  ";cell id; ts(new|bhabha) - ts(previous iteration|merged) [ns]", 8736,
85  1, 8736 + 1));
86  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld__crateID(new TH1F("tcrateNew_MINUS_tcrateOld__crateID",
87  ";crate id; tcrate(new | bhabha) - tcrate(previous iteration | merged) [ns]",
88  52, 1, 52 + 1));
89  unique_ptr<TH1F> tsNew_MINUS_tsCustomPrev__cid(new TH1F("TsNew_MINUS_TsCustomPrev__cid",
90  ";cell id; ts(new|bhabha) - ts(old = 'before 1st iter'|merged) [ns]",
91  8736, 1, 8736 + 1));
92  unique_ptr<TH1F> tsNew_MINUS_tsOldBhabha__cid(new TH1F("TsNew_MINUS_TsOldBhabha__cid",
93  ";cell id; ts(new|bhabha) - ts(previous iteration|bhabha) [ns]", 8736,
94  1, 8736 + 1));
95 
96 
97  // Histogram of the new time constant values minus values from previous iteration
98  unique_ptr<TH1F> tsNew_MINUS_tsOld(new TH1F("TsNew_MINUS_TsOld",
99  ";ts(new | bhabha) - ts(previous iteration | merged) [ns];Number of crystals",
100  201, -10.05, 10.05));
101  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld(new TH1F("tcrateNew_MINUS_tcrateOld",
102  ";tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
103  201, -10.05, 10.05));
104  unique_ptr<TH1F> tsNew_MINUS_tsCustomPrev(new TH1F("TsNew_MINUS_TsCustomPrev",
105  ";ts(new | bhabha) - ts(old = 'before 1st iter' | merged) [ns];Number of crystals",
106  285, -69.5801, 69.5801));
107  unique_ptr<TH1F> tsNew_MINUS_tsOldBhabha(new TH1F("TsNew_MINUS_TsOldBhabha",
108  ";ts(new | bhabha) - ts(previous iteration | bhabha) [ns];Number of crystals",
109  201, -10.05, 10.05));
110 
111 
112 
113  // Histograms just for crates and collecting for all the different runs rather than run by run.
114  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld_allRuns(new TH1F("tcrateNew_MINUS_tcrateOld_allRuns",
115  "Crate time constant changes over all runs : tcrate(new) uncertainty < 0.1ns;tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
116  201, -10.05, 10.05));
117 
118  unique_ptr<TH1F> tcrateNew_MINUS_tcrateOld_allRuns_allCrates(new TH1F("tcrateNew_MINUS_tcrateOld_allRuns_allCrates",
119  "Crate time constant changes over all runs : all crates;tcrate(new) - tcrate(previous iteration) [ns];Number of crates",
120  201, -10.05, 10.05));
121 
122  unique_ptr<TH1I> num_tcrates_perRun(new TH1I("num_tcrates_perRun",
123  "Number of good tcrates in each run;Run number;Number of good tcrates",
124  6000, 0, 6000));
125 
126  unique_ptr<TH2F> tcrateNew_MINUS_tcrateOld__vs__runNum(new TH2F("tcrateNew_MINUS_tcrateOld__vs__runNum",
127  "Crate time constant changes vs run number : tcrate(new) uncertainty < 0.1ns;Run number;tcrate(new) - tcrate(previous iteration) [ns]",
128  6000, 0, 6000, 21, -10.5, 10.5));
129 
130 
131 
132 
133 
134  if (!TimevsCrysNoCalibrations) return c_Failure;
135 
138  TFile* histfile = 0;
139  TFile* histExtraCrateInfofile = 0; // Only created if needed for crates
140 
141 
142  unique_ptr<TTree> tree_crystal(new TTree("tree_crystal", "Debug data from bhabha time calibration algorithm for crystals"));
143 
144  unique_ptr<TTree> tree_crate(new TTree("tree_crate", "Debug data from bhabha time calibration algorithm for crates"));
145  int tree_cid;
146 
147  // Vector of time offsets to be saved in the database.
148  vector<float> t_offsets;
149  // Vector of time offset uncertainties to be saved in the database.
150  vector<float> t_offsets_unc;
151 
152  vector<float> t_offsets_prev; // previous time offsets
153 
154  // temp copies of offsets
155  //vector<float> t_offsets_temp;
156  //vector<float> t_offsets_unc_temp;
157 
158 
159  int minNumEntries = 40;
160  int minNumEntriesCrateConvergence = 1000;
161 
162 
163  double mean = 0;
164  double sigma = -1;
165  double mean_unc = 0;
166  int crystalCalibSaved = 0;
167  double tsPrev = 0;
168 
169 
170  bool minRunNumBool = false;
171  bool maxRunNumBool = false;
172  int minRunNum = -1;
173  int maxRunNum = -1;
174  int minExpNum = -1;
175  int maxExpNum = -1;
176  for (auto expRun : getRunList()) {
177  int expNumber = expRun.first;
178  int runNumber = expRun.second;
179  if (!minRunNumBool) {
180  minExpNum = expNumber;
181  minRunNum = runNumber;
182  minRunNumBool = true;
183  }
184  if (!maxRunNumBool) {
185  maxExpNum = expNumber;
186  maxRunNum = runNumber;
187  maxRunNumBool = true;
188  }
189  if (((minRunNum > runNumber) && (minExpNum >= expNumber)) ||
190  (minExpNum > expNumber)) {
191  minExpNum = expNumber;
192  minRunNum = runNumber;
193  }
194  if (((maxRunNum < runNumber) && (maxExpNum <= expNumber)) ||
195  (maxExpNum < expNumber))
196 
197  {
198  maxExpNum = expNumber;
199  maxRunNum = runNumber;
200  }
201  }
202 
203  B2INFO("debugFilenameBase = " << debugFilenameBase);
204  string runNumsString = string("_") + to_string(minExpNum) + "_" + to_string(minRunNum) + string("-") +
205  to_string(maxExpNum) + "_" + to_string(maxRunNum);
206  string debugFilename = debugFilenameBase + runNumsString + string(".root");
207  string extraCratedebugFilename = debugFilenameBase + string("_cratesAllRuns.root");
208 
209 
210  // Need to load information about the event/run/experiment to get the right database information
211  // Will be used for:
212  // * ECLChannelMapper (to map crystal to crates)
213  // * crystal payload updating for iterating crystal and crate fits
214  int eventNumberForCrates = 1;
215 
217  // simulate the initialize() phase where we can register objects in the DataStore
219  evtPtr.registerInDataStore();
221  // now construct the event metadata
222  evtPtr.construct(eventNumberForCrates, minRunNum, minExpNum);
223  // and update the database contents
224  DBStore& dbstore = DBStore::Instance();
225  dbstore.update();
226  // this is only needed it the payload might be intra-run dependent,
227  // that is if it might change during one run as well
228  dbstore.updateEvent();
229 
230 
231  B2INFO("Uploading payload for exp " << minExpNum << ", run " << minRunNum << ", event " << eventNumberForCrates);
232  updateDBObjPtrs(eventNumberForCrates, minRunNum, minExpNum);
233  unique_ptr<ECLChannelMapper> crystalMapper(new ECL::ECLChannelMapper());
234  crystalMapper->initFromDB();
235 
236 
237  //------------------------------------------------------------------------
238  //..Read payloads from database
239  B2INFO("Reading payloads: ECLCrystalTimeOffset and ECLCrateTimeOffset");
240  DBObjPtr<Belle2::ECLCrystalCalib> crystalTimeObject("ECLCrystalTimeOffset");
241  DBObjPtr<Belle2::ECLCrystalCalib> crateTimeObject("ECLCrateTimeOffset");
242 
243  //..Get vectors of values from the payloads
244  vector<float> currentValuesCrys = crystalTimeObject->getCalibVector();
245  vector<float> currentUncCrys = crystalTimeObject->getCalibUncVector();
246  vector<float> currentValuesCrate = crateTimeObject->getCalibVector();
247  vector<float> currentUncCrate = crateTimeObject->getCalibUncVector();
248 
249  //..Print out a few values for quality control
250  B2INFO("Values read from database. Write out for their values for comparison against those from tcol");
251  for (int ic = 0; ic < 8736; ic += 500) {
252  B2INFO("ts: cellID " << ic + 1 << " " << currentValuesCrys[ic] << " +/- " << currentUncCrys[ic]);
253  B2INFO("tcrate: cellID " << ic + 1 << " " << currentValuesCrate[ic] << " +/- " << currentUncCrate[ic]);
254  }
255 
256 
257  //..Read in the previous crystal payload values
258  vector<float> prevValuesCrys(8736);
259  if (readPrevCrysPayload) {
260  DBObjPtr<Belle2::ECLCrystalCalib> customPrevCrystalTimeObject("ECLCrystalTimeOffsetPreviousValues");
261  //..Get vectors of values from the payloads
262  prevValuesCrys = customPrevCrystalTimeObject->getCalibVector();
263 
264  //..Print out a few values for quality control
265  B2INFO("Previous values read from database. Write out for their values for comparison against those from tcol");
266  for (int ic = 0; ic < 8736; ic += 500) {
267  B2INFO("ts custom previous payload: cellID " << ic + 1 << " " << prevValuesCrys[ic]);
268  }
269  }
270 
271 
272  //..Read bhabha payloads from database
273  B2INFO("Reading payloads: ECLCrystalTimeOffsetBhabha");
274  DBObjPtr<Belle2::ECLCrystalCalib> crystalBhabhaTimeObject("ECLCrystalTimeOffsetBhabha");
275 
276  //..Get vectors of values from the payloads
277  vector<float> currentBhabhaValuesCrys = crystalBhabhaTimeObject->getCalibVector();
278  vector<float> currentBhabhaUncCrys = crystalBhabhaTimeObject->getCalibUncVector();
279 
280  //..Print out a few values for quality control
281  for (int ic = 0; ic < 8736; ic += 500) {
282  B2INFO("ts bhabha: cellID " << ic + 1 << " " << currentBhabhaValuesCrys[ic] << " +/- " << currentBhabhaUncCrys[ic]);
283  }
284 
285 
286  //------------------------------------------------------------------------
287  //..Get the reference crystal information
288  auto refCrysIDzeroingCrate = getObjectPtr<TH1F>("refCrysIDzeroingCrate");
289 
290  // crystal index for the crystals (one per crate) that is used as the reference crystal. This one has
291  // ts defined as zero. The crystal id runs 1...8636, not starting at 0.
292  B2INFO("Extract reference crystals from collector histogram.");
293  vector <short> crystalIDreferenceUntested;
294  for (int bin = 1; bin <= 8736; bin++) {
295  if (refCrysIDzeroingCrate->GetBinContent(bin) > 0.5) {
296  crystalIDreferenceUntested.push_back(bin);
297  }
298  }
299 
300  // Output for the user the crystal id to be used as a reference
301  // and also state which crate the crystal is in.
302  B2INFO("Reference crystals to define as having ts=0. Base 1 counting for both crates and crystals");
303  for (long unsigned int crysRefCounter = 0; crysRefCounter < crystalIDreferenceUntested.size(); crysRefCounter++) {
304  int crys_id = crystalIDreferenceUntested[crysRefCounter] ;
305  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
306  B2INFO(" crystal " << crys_id << " is a reference for crate " << crate_id_from_crystal);
307  }
308 
309  // Check that the reference crystals make sense. There should be exactly one per crate but
310  // since the crystal calibration executes over multiple runs, it is possible that the
311  // reference crystals could change but we don't want to allow this.
312  B2INFO("Checking number of reference crystals");
313  B2INFO("Number of reference crystals = " << crystalIDreferenceUntested.size());
314 
315  // Check that there are exactly 52 reference crystals
316  if (crystalIDreferenceUntested.size() != 52) {
317  B2FATAL("The number of reference crystals does not equal 52, which is one per crate");
318  return c_Failure;
319  } else {
320  B2INFO("Number of reference crystals is 52 as required");
321  }
322 
323  // Count the number of reference crystals for each crate to make sure that there is exactly
324  // one reference crystal for each crate as defined by the payload/database.
325  // also fill the final vector that maps the crate id to the reference crystal id
326  vector <short> crateIDsNumRefCrystalsUntested(52, 0);
327  vector <short> crystalIDReferenceForZeroTs(52, 0);
328 
329  for (long unsigned int crysRefCounter = 0; crysRefCounter < crystalIDreferenceUntested.size(); crysRefCounter++) {
330  int crys_id = crystalIDreferenceUntested[crysRefCounter] ;
331  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
332  crateIDsNumRefCrystalsUntested[crate_id_from_crystal - 1]++;
333  crystalIDReferenceForZeroTs[crate_id_from_crystal - 1] = crys_id;
334  }
335  B2INFO("crystalIDReferenceForZeroTs = {");
336  for (int crateTest = 0; crateTest < 52 - 1; crateTest++) {
337  B2INFO(crystalIDReferenceForZeroTs[crateTest] << ",");
338  }
339  B2INFO(crystalIDReferenceForZeroTs[52 - 1] << "}");
340 
341 
342  // Make sure that there is only one reference crystal per crate as defined by the payload/database
343  for (int crateTest = 0; crateTest < 52; crateTest++) {
344  if (crateIDsNumRefCrystalsUntested[crateTest] != 1) {
345  B2FATAL("Crate " << crateTest + 1 << " (base 1) has " << crateIDsNumRefCrystalsUntested[crateTest] << " reference crystals");
346  return c_Failure;
347  }
348  }
349  B2INFO("All reference crystals are reasonably mapped one crystal to one crate for all crates");
350 
351 
352 
353  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.");
354 
355  /* Test if the user wants to modify the reference crystals. This will probably be done very rarely,
356  perhaps less than once per year. If the user does want to change one or more reference
357  crystals then perform the checks again to make sure that there is still just one reference crystal
358  per crate after the modifications to the payload values with the user values.*/
359  bool userSetRefCrysPerCrate = false ;
360  for (int crateTest = 0; crateTest < 52; crateTest++) {
361  if (refCrysPerCrate[crateTest] != -1) {
362  crystalIDReferenceForZeroTs[crateTest] = refCrysPerCrate[crateTest] ;
363  B2INFO("Crate " << crateTest + 1 << " (base 1) new reference crystal = " << crystalIDReferenceForZeroTs[crateTest]);
364  userSetRefCrysPerCrate = true ;
365  }
366  }
367  if (userSetRefCrysPerCrate) {
368  B2INFO("User changed reference crystals via steering script");
369 
370  // Validate crystals per crate again with the new user set values
371  fill(crateIDsNumRefCrystalsUntested.begin(), crateIDsNumRefCrystalsUntested.end(), 0);
372  for (long unsigned int crysRefCounter = 0; crysRefCounter < 52; crysRefCounter++) {
373  int crys_id = crystalIDReferenceForZeroTs[crysRefCounter] ;
374  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
375  crateIDsNumRefCrystalsUntested[crate_id_from_crystal - 1]++;
376  }
377  for (int crateTest = 0; crateTest < 52; crateTest++) {
378  if (crateIDsNumRefCrystalsUntested[crateTest] != 1) {
379  B2FATAL("Crate " << crateTest + 1 << " (base 1) has " << crateIDsNumRefCrystalsUntested[crateTest] << " reference crystals");
380  return c_Failure;
381  }
382  }
383  B2INFO("All reference crystals are reasonably mapped one crystal to one crate for all crates after changes made by user steering script.");
384 
385  // Save the information to the payload if there is at least one crate
386  // reference crystal that has been modified by the user steering file
388  refCrystalsCalib->setCalibVector(crystalIDReferenceForZeroTs);
389  saveCalibration(refCrystalsCalib, "ECLReferenceCrystalPerCrateCalib");
390  B2INFO("Created reference crystal per crate payload: ECLReferenceCrystalPerCrateCalib");
391  } else {
392  B2INFO("User did not change reference crystals via steering script");
393  }
394 
395 
396  //------------------------------------------------------------------------
397  //..Start looking at timing information
398 
399  B2INFO("Debug output rootfile: " << debugFilename);
400  histfile = new TFile(debugFilename.c_str(), "recreate");
401 
402 
403  TimevsCrysPrevCrateCalibPrevCrystCalib ->Write();
404  TimevsCratePrevCrateCalibPrevCrystCalib->Write();
405  TimevsCrysNoCalibrations ->Write();
406  TimevsCrysPrevCrateCalibNoCrystCalib ->Write();
407  TimevsCrateNoCrateCalibPrevCrystCalib ->Write();
408 
409  cutflow->Write();
410 
411 
412  if (debugOutput) {
413  tree_crystal->Branch("cid", &tree_cid)->SetTitle("Cell ID, 1..8736");
414  tree_crystal->Branch("ts", &mean)->SetTitle("Time offset mean, ts, ns");
415  tree_crystal->Branch("tsUnc", &mean_unc)->SetTitle("Error of time ts mean, ns.");
416  tree_crystal->Branch("tsSigma", &sigma)->SetTitle("Sigma of time ts distribution, ns");
417  tree_crystal->Branch("crystalCalibSaved",
418  &crystalCalibSaved)->SetTitle("0=crystal skipped, 1=crystal calib saved (num entries based)");
419  tree_crystal->Branch("tsPrev", &tsPrev)->SetTitle("Previous crystal time offset, ts, ns");
420  tree_crystal->SetAutoSave(10);
421  }
422 
423 
424  double hist_tmin = TimevsCrysNoCalibrations->GetYaxis()->GetXmin();
425  double hist_tmax = TimevsCrysNoCalibrations->GetYaxis()->GetXmax();
426 
427  double time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
428  double time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
429 
430  B2INFO("hist_tmin = " << hist_tmin);
431  B2INFO("hist_tmax = " << hist_tmax);
432 
433  /* 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency, fRF=508.889 MHz.
434  Same for all crystals. Proper accurate value*/
435  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::m_rf) * 1e3;
436 
437  // The ts and tcrate database values are filled once per tcol instance so count the number of times that the database values
438  // were summed together by the histogram merging process and extract out the original values again.
439  auto databaseCounter = getObjectPtr<TH1I>("databaseCounter");
440  float numTimesFilled = databaseCounter->GetBinContent(1);
441  B2INFO("Number of times database histograms were merged = " << numTimesFilled);
442 
443 
444  auto TsDatabase = getObjectPtr<TH1F>("TsDatabase");
445  auto TsDatabaseUnc = getObjectPtr<TH1F>("TsDatabaseUnc");
446  for (int i = 1; i <= 8736; i++) {
447  t_offsets.push_back(TsDatabase->GetBinContent(i) / numTimesFilled);
448  t_offsets_prev.push_back(TsDatabase->GetBinContent(i) / numTimesFilled);
449 
450  B2INFO("t_offsets_prev (last iter) at crysID " << i << " = " << t_offsets_prev[i - 1]);
451 
452  t_offsets_unc.push_back(TsDatabaseUnc->GetBinContent(i) / numTimesFilled);
453  }
454 
455 
456  /* CRYSTAL CORRECTIONS */
457 
458  /* Make a 1D histogram of the number of hits per crystal. This will help with the validations
459  to make sure that all the crystals had enough hits and to look for problems.*/
460  TH1D* h_crysHits = TimevsCrysPrevCrateCalibNoCrystCalib->ProjectionX("h_crysHits");
461  h_crysHits->SetTitle("Hits per crystal;Crystal id");
462 
463  histfile->WriteTObject(h_crysHits, "h_crysHits");
464 
465 
466  // Loop over all the crystals for doing the crystal calibation
467  for (int crys_id = cellIDLo; crys_id <= cellIDHi; crys_id++) {
468  mean = 0;
469  mean_unc = -1;
470  sigma = -1;
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  double database_mean_crate = 0;
804  double database_mean_crate_unc = 0;
805 
806  double mean_crate = 0;
807  double sigma_crate = -1;
808 
809  vector<float> tcrate_mean_new(52, 0.0);
810  vector<float> tcrate_mean_unc_new(52, 0.0);
811  vector<float> tcrate_sigma_new(52, 0.0);
812  vector<float> tcrate_mean_prev(52, 0.0);
813  vector<float> tcrate_sigma_prev(52, 0.0);
814  vector<bool> tcrate_new_was_set(52, false);
815  vector<bool> tcrate_new_goodQuality(52, false);
816 
817  B2DEBUG(22, "crate vectors set");
818 
819 
820 
821  // Crate time calibration constants are saved per crystal so read them per crystal
822  // and save as one entry per crate in the array
823  for (int crys_id = 1; crys_id <= 8736; crys_id++) {
824  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
825  tcrate_mean_prev[crate_id_from_crystal - 1] = TcrateDatabase->GetBinContent(crys_id) / numTimesFilled;
826  }
827 
828 
829  B2INFO("Print out previous crate time calibration constants to make sure they match from the two different sources.");
830  for (int crate_id = 1; crate_id <= 52; crate_id++) {
831  B2INFO("tcrate_mean_prev[crate " << crate_id << " (base 1)] = " << tcrate_mean_prev[crate_id - 1]);
832 
833  int thisRefCellID = crystalIDReferenceForZeroTs[crate_id - 1];
834  B2INFO("tcrate from payload: ref cellID " << thisRefCellID << " " << currentValuesCrate[thisRefCellID - 1] << " +/- " <<
835  currentUncCrate[thisRefCellID - 1]);
836  }
837 
838 
839  /* Read in the histogram about the crate calibration constant differences between iterations
840  if it exists. This way the histograms can be updated after each run.*/
841  if (crateIDLo <= crateIDHi) {
842  TFile* histExtraCrateInfofile_dummy = 0; // Only created if needed for crates
843  B2INFO("Debug output rootfile used for crate iterations: " << extraCratedebugFilename);
844  histExtraCrateInfofile_dummy = new TFile(extraCratedebugFilename.c_str(), "UPDATE");
845 
846 
847  // If the histogram already exists then read it in and recreate the file so that we can save an updated version of the histogram.
848  TKey* key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld_allRuns");
849  if (key != 0) {
850  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld_allRuns");
851  tcrateNew_MINUS_tcrateOld_allRuns->Add(h);
852  }
853 
854  key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
855  if (key != 0) {
856  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
857  tcrateNew_MINUS_tcrateOld_allRuns_allCrates->Add(h);
858  }
859 
860  key = histExtraCrateInfofile_dummy->FindKey("num_tcrates_perRun");
861  if (key != 0) {
862  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("num_tcrates_perRun");
863  num_tcrates_perRun->Add(h);
864  }
865 
866  key = histExtraCrateInfofile_dummy->FindKey("tcrateNew_MINUS_tcrateOld__vs__runNum");
867  if (key != 0) {
868  TH1F* h = (TH1F*)histExtraCrateInfofile_dummy->Get("tcrateNew_MINUS_tcrateOld__vs__runNum");
869  tcrateNew_MINUS_tcrateOld__vs__runNum->Add(h);
870  }
871 
872  histExtraCrateInfofile_dummy->Close();
873 
874  /* After reading in all the histograms, recreate the root file from empty so that the
875  histograms can be made again with the updated values.*/
876  histExtraCrateInfofile = new TFile(extraCratedebugFilename.c_str(), "recreate");
877  }
878 
879 
880 
881  for (int crate_id = crateIDLo; crate_id <= crateIDHi; crate_id++) {
882 
883  B2DEBUG(22, "Start of crate id = " << crate_id);
884 
885  TH1D* h_time_crate = TimevsCrateNoCrateCalibPrevCrystCalib->ProjectionY("h_time_psi_crate", crate_id, crate_id);
886  TH1D* h_time_crate_mask = (TH1D*)h_time_crate->Clone();
887  TH1D* h_time_crate_masked = (TH1D*)h_time_crate->Clone();
888  TH1D* h_time_crate_rebin = (TH1D*)h_time_crate->Clone();
889 
890 
891  // Do rebinning and cleaning of some bins but only if the user selection values call for it since it slows the code down
892  if (meanCleanRebinFactor != 1 || meanCleanCutMinFactor != 1) {
893 
894  h_time_crate_rebin->Rebin(meanCleanRebinFactor);
895  h_time_crate_mask->Scale(0.0); // set all bins to being masked off
896 
897  time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
898  time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
899 
900  // Find value of bin with max value
901  double histRebin_max = h_time_crate_rebin->GetMaximum();
902 
903  bool maskedOutNonZeroBin = false;
904  // Loop over all bins to find those with content less than a certain threshold. Mask the non-rebinned histogram for the corresponding bins
905  for (int bin = 1; bin <= h_time_crate_rebin->GetNbinsX(); bin++) {
906  for (int rebinCounter = 1; rebinCounter <= meanCleanRebinFactor; rebinCounter++) {
907  int nonRebinnedBinNumber = (bin - 1) * meanCleanRebinFactor + rebinCounter;
908  if (nonRebinnedBinNumber < h_time_crate->GetNbinsX()) {
909  if (h_time_crate_rebin->GetBinContent(bin) >= histRebin_max * meanCleanCutMinFactor) {
910  h_time_crate_mask->SetBinContent(nonRebinnedBinNumber, 1);
911 
912  // Save the lower and upper edges of the rebin histogram time range for fitting purposes
913  double x_lower = h_time_crate_rebin->GetXaxis()->GetBinLowEdge(bin);
914  double x_upper = h_time_crate_rebin->GetXaxis()->GetBinUpEdge(bin);
915  if (x_lower < time_fit_min) {
916  time_fit_min = x_lower;
917  }
918  if (x_upper > time_fit_max) {
919  time_fit_max = x_upper;
920  }
921  } else {
922  if (h_time_crate->GetBinContent(nonRebinnedBinNumber) > 0) {
923  B2DEBUG(22, "Setting bin " << nonRebinnedBinNumber << " from " << h_time_crate_masked->GetBinContent(
924  nonRebinnedBinNumber) << " to 0");
925  maskedOutNonZeroBin = true;
926  }
927  h_time_crate_masked->SetBinContent(nonRebinnedBinNumber, 0);
928  }
929  }
930  }
931  }
932  B2INFO("Bins with non-zero values have been masked out: " << maskedOutNonZeroBin);
933  h_time_crate_masked->ResetStats();
934  h_time_crate_mask->ResetStats();
935 
936  }
937 
938 
939 
940  B2DEBUG(22, "crate loop - projected h_time_psi_crate");
941 
942 
943  double default_mean_crate = h_time_crate_masked->GetMean();
944  double default_mean_crate_unc = h_time_crate_masked->GetMeanError();
945  double default_sigma_crate = h_time_crate_masked->GetStdDev();
946  B2INFO("Fitting crate between " << time_fit_min << " and " << time_fit_max);
947  TF1* gaus = new TF1("func", "gaus(0)", time_fit_min, time_fit_max);
948  gaus->SetParNames("numCrateHisNormalization", "mean", "sigma");
949  double hist_max = h_time_crate->GetMaximum();
950  double stddev = h_time_crate->GetStdDev();
951  sigma_crate = stddev;
952  mean_crate = default_mean_crate;
953  gaus->SetParameter(0, hist_max / 2.);
954  gaus->SetParameter(1, mean_crate);
955  gaus->SetParameter(2, sigma_crate);
956 
957  h_time_crate_masked->Fit(gaus, "LQR"); // L for likelihood, R for x-range, Q for fit quiet mode
958 
959  double fit_mean_crate = gaus->GetParameter(1);
960  double fit_mean_crate_unc = gaus->GetParError(1);
961  double fit_sigma_crate = gaus->GetParameter(2);
962 
963  double meanDiff = fit_mean_crate - default_mean_crate;
964  double meanUncDiff = fit_mean_crate_unc - default_mean_crate_unc;
965  double sigmaDiff = fit_sigma_crate - default_sigma_crate;
966 
967  bool good_fit = false;
968 
969  B2DEBUG(22, "Crate id = " << crate_id << " with crate mean = " << default_mean_crate << " +- " << fit_mean_crate_unc);
970 
971  if ((fabs(meanDiff) > 7) ||
972  (fabs(meanUncDiff) > 7) ||
973  (fabs(sigmaDiff) > 7) ||
974  (fit_mean_crate_unc > 3) ||
975  (fit_sigma_crate < 0.1) ||
976  (fit_mean_crate < time_fit_min) ||
977  (fit_mean_crate > time_fit_max)) {
978  B2INFO("Crate id = " << crate_id);
979  B2INFO("fit mean, default mean = " << fit_mean_crate << ", " << default_mean_crate);
980  B2INFO("fit sigma, default sigma = " << fit_sigma_crate << ", " << default_sigma_crate);
981 
982  B2INFO("crate fit mean - hist mean = " << meanDiff);
983  B2INFO("fit mean unc. - hist mean unc. = " << meanUncDiff);
984  B2INFO("fit sigma - hist sigma = " << sigmaDiff);
985  B2INFO("fit_mean_crate = " << fit_mean_crate);
986  B2INFO("time_fit_min = " << time_fit_min);
987  B2INFO("time_fit_max = " << time_fit_max);
988  } else {
989  good_fit = true;
990  }
991 
992  int numEntries = h_time_crate->GetEntries();
993  B2INFO("Number entries = " << numEntries);
994  if ((numEntries >= minNumEntries) && good_fit) {
995  database_mean_crate = fit_mean_crate;
996  database_mean_crate_unc = fit_mean_crate_unc;
997 
998  if ((numEntries >= minNumEntriesCrateConvergence) && (fit_mean_crate_unc < 0.1)) {
999  tcrate_new_goodQuality[crate_id - 1] = true;
1000  }
1001  } else {
1002  database_mean_crate = default_mean_crate;
1003  database_mean_crate_unc = default_mean_crate_unc;
1004  }
1005  if (numEntries == 0) {
1006  database_mean_crate = 0;
1007  database_mean_crate_unc = 0;
1008  }
1009 
1010  tcrate_mean_new[crate_id - 1] = database_mean_crate;
1011  tcrate_mean_unc_new[crate_id - 1] = database_mean_crate_unc;
1012  tcrate_sigma_new[crate_id - 1] = fit_sigma_crate;
1013  tcrate_new_was_set[crate_id - 1] = true;
1014 
1015 
1016  histfile->WriteTObject(h_time_crate, (string("h_time_psi_crate") + to_string(crate_id)).c_str());
1017  histfile->WriteTObject(h_time_crate_masked, (string("h_time_psi_crate_masked") + to_string(crate_id)).c_str());
1018  histfile->WriteTObject(h_time_crate_rebin, (string("h_time_psi_crate_rebinned") + to_string(crate_id)).c_str());
1019 
1020  delete gaus;
1021  }
1022 
1023  B2DEBUG(22, "crate histograms made");
1024 
1025 
1026  // Save database for crates
1027  // Vector of time offsets to be saved in the database.
1028  vector<float> t_offsets_crate;
1029  // Vector of time offset uncertainties to be saved in the database.
1030  vector<float> t_offsets_crate_unc;
1031  for (int i = 1; i <= 8736; i++) {
1032  t_offsets_crate.push_back(0);
1033  t_offsets_crate_unc.push_back(0);
1034  }
1035 
1036 
1037  for (int crys_id = 1; crys_id <= 8736; crys_id++) {
1038  int crate_id_from_crystal = crystalMapper->getCrateID(crys_id);
1039  if (tcrate_new_was_set[crate_id_from_crystal - 1]) {
1040  t_offsets_crate[crys_id - 1] = tcrate_mean_new[crate_id_from_crystal - 1] / TICKS_TO_NS;
1041  t_offsets_crate_unc[crys_id - 1] = tcrate_mean_unc_new[crate_id_from_crystal - 1] / TICKS_TO_NS;
1042 
1043  } else {
1044  t_offsets_crate[crys_id - 1] = tcrate_mean_prev[crate_id_from_crystal - 1];
1045  B2INFO("used old crate mean but zeroed uncertainty since not saved in root file");
1046  }
1047  }
1048 
1049 
1050  // Fill histograms with the difference in the tcrate values
1051  if (crateIDLo <= crateIDHi) {
1052  for (int crate_id = crateIDLo; crate_id <= crateIDHi; crate_id++) {
1053  // tcrate_mean_new already in ns, but tcrate_mean_prev in ticks.
1054  double tCrateDiff_ns = tcrate_mean_new[crate_id - 1] - (tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS);
1055  B2INFO("Crate " << crate_id << ": tcrate new - previous iteration = "
1056  << tcrate_mean_new[crate_id - 1]
1057  << " - " << tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS
1058  << " = " << tCrateDiff_ns << " ns");
1059  tcrateNew_MINUS_tcrateOld__crateID->SetBinContent(crate_id, tCrateDiff_ns);
1060  tcrateNew_MINUS_tcrateOld__crateID->SetBinError(crate_id, 0);
1061  tcrateNew_MINUS_tcrateOld__crateID->ResetStats();
1062 
1063  tcrateNew_MINUS_tcrateOld->Fill(tCrateDiff_ns);
1064 
1065  // Save the histograms monitoring the change between iterations
1066  tcrateNew_MINUS_tcrateOld_allRuns_allCrates->Fill(tCrateDiff_ns);
1067  if (tcrate_new_goodQuality[crate_id - 1]) {
1068  tcrateNew_MINUS_tcrateOld_allRuns->Fill(tCrateDiff_ns);
1069  num_tcrates_perRun->Fill(minRunNum);
1070  tcrateNew_MINUS_tcrateOld__vs__runNum->Fill(minRunNum, tCrateDiff_ns);
1071  }
1072  }
1073 
1074  // Save the histograms to the output root file
1075  histfile->WriteTObject(tcrateNew_MINUS_tcrateOld__crateID.get(), "tcrateNew_MINUS_tcrateOld__crateID");
1076  histfile->WriteTObject(tcrateNew_MINUS_tcrateOld.get(), "tcrateNew_MINUS_tcrateOld");
1077 
1078  // Save the histograms ot the crate iterations file
1079  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld_allRuns.get(), "tcrateNew_MINUS_tcrateOld_allRuns");
1080  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld_allRuns_allCrates.get(),
1081  "tcrateNew_MINUS_tcrateOld_allRuns_allCrates");
1082  histExtraCrateInfofile->WriteTObject(num_tcrates_perRun.get(), "num_tcrates_perRun");
1083  histExtraCrateInfofile->WriteTObject(tcrateNew_MINUS_tcrateOld__vs__runNum.get(), "tcrateNew_MINUS_tcrateOld__vs__runNum");
1084  }
1085 
1086 
1087  ECLCrystalCalib* BhabhaTCrateCalib = new ECLCrystalCalib();
1088  BhabhaTCrateCalib->setCalibVector(t_offsets_crate, t_offsets_crate_unc);
1089 
1090 
1091  // Save the information to the payload if there is at least one crate
1092  // begin calibrated.
1093  if (crateIDLo <= crateIDHi) {
1094  saveCalibration(BhabhaTCrateCalib, "ECLCrateTimeOffset");
1095  B2DEBUG(22, "crate payload made");
1096 
1097  histExtraCrateInfofile->Close();
1098  }
1099 
1100 
1101  int tree_crateid;
1102  int tree_runNum;
1103  double tree_tcrate_mean;
1104  double tree_tcrate_mean_unc;
1105  double tree_tcrate_sigma;
1106  double tree_tcrate_meanPrev;
1107 
1108  tree_crate->Branch("runNum", &tree_runNum)->SetTitle("Run number, 0..infinity and beyond!");
1109  tree_crate->Branch("crateid", &tree_crateid)->SetTitle("Crate id, 1..52");
1110  tree_crate->Branch("tcrate", &tree_tcrate_mean)->SetTitle("Crate time offset mean, tcrate, ns");
1111  tree_crate->Branch("tcratePrev", &tree_tcrate_meanPrev)->SetTitle("Previous crate time offset mean, tcrate, ns");
1112  tree_crate->Branch("tcrate_unc", &tree_tcrate_mean_unc)->SetTitle("Error of time tcrate mean, ns.");
1113  tree_crate->Branch("tcrate_sigma", &tree_tcrate_sigma)->SetTitle("Sigma of time tcrate distribution, ns");
1114  tree_crate->SetAutoSave(10);
1115 
1116 
1117  for (auto expRun : getRunList()) {
1118  // Key command to make sure your DBObjPtrs are correct
1119  B2INFO("run num, exp num: " << expRun.second << ", " << expRun.first);
1120  int runNumber = expRun.second;
1121 
1122  for (int crate_id = 1; crate_id <= 52; crate_id++) {
1123  if (tcrate_new_was_set[crate_id - 1]) {
1124  tree_runNum = runNumber;
1125  tree_crateid = crate_id;
1126  tree_tcrate_mean = tcrate_mean_new[crate_id - 1];
1127  tree_tcrate_mean_unc = tcrate_mean_unc_new[crate_id - 1];
1128  tree_tcrate_sigma = tcrate_sigma_new[crate_id - 1];
1129  tree_tcrate_meanPrev = tcrate_mean_prev[crate_id - 1] * TICKS_TO_NS;
1130  tree_crate->Fill();
1131  }
1132  }
1133  }
1134 
1135  B2DEBUG(22, "end of crate corrections .....");
1136 
1137  tree_crystal->Write();
1138  tree_crate->Write();
1139 
1140  histfile->Close();
1141 
1142  tree_crystal = 0;
1143  tree_crate = 0;
1144 
1145  B2INFO("Finished talgorithm");
1146  return c_OK;
1147 }
1148 
Belle2::ECL::eclBhabhaTAlgorithm::readPrevCrysPayload
bool readPrevCrysPayload
Read the previous crystal payload values for comparison.
Definition: eclBhabhaTAlgorithm.h:59
Belle2::ECL::eclBhabhaTAlgorithm::refCrysPerCrate
int refCrysPerCrate[52]
List of crystals, one per crate, used as reference time for crystal time calibration.
Definition: eclBhabhaTAlgorithm.h:64
Belle2::DataStore::Instance
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
Belle2::ECL::eclBhabhaTAlgorithm::cellIDHi
int cellIDHi
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
Definition: eclBhabhaTAlgorithm.h:51
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::ECL::ECLChannelMapper::getCrateID
int getCrateID(int iCOPPERNode, int iFINESSE)
get crate number by given COPPER node number and FINESSE number
Definition: ECLChannelMapper.cc:211
Belle2::ECL::eclBhabhaTAlgorithm::debugOutput
bool debugOutput
Save every histogram and fitted function to debugFilename.
Definition: eclBhabhaTAlgorithm.h:60
Belle2::DataStore::setInitializeActive
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
Belle2::ECLReferenceCrystalPerCrateCalib::setCalibVector
void setCalibVector(const std::vector< short > &refCrystals)
Set vector of constants with uncertainties.
Definition: ECLReferenceCrystalPerCrateCalib.h:62
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::ECLCrystalCalib
General DB object to store one calibration number per ECL crystal.
Definition: ECLCrystalCalib.h:34
Belle2::ECL::eclBhabhaTAlgorithm::crateIDLo
int crateIDLo
Fit crates with crateID0 in the inclusive range [crateIDLo,crateIDHi].
Definition: eclBhabhaTAlgorithm.h:56
Belle2::ECL::eclBhabhaTAlgorithm::savePrevCrysPayload
bool savePrevCrysPayload
Save the previous crystal payload values for comparison.
Definition: eclBhabhaTAlgorithm.h:58
Belle2::ECL::eclBhabhaTAlgorithm::calibrate
EResult calibrate() override
..Run algorithm on events
Definition: eclBhabhaTAlgorithm.cc:46
Belle2::CalibrationAlgorithm::updateDBObjPtrs
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
Definition: CalibrationAlgorithm.cc:398
Belle2::DBObjPtr< Belle2::ECLCrystalCalib >
Belle2::ECL::eclBhabhaTAlgorithm::crateIDHi
int crateIDHi
Fit crates with crateID0 in the inclusive range [crateIDLo,crateIDHi].
Definition: eclBhabhaTAlgorithm.h:57
Belle2::ECL::ECLChannelMapper
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
Definition: ECLChannelMapper.h:36
Belle2::DBStore
Singleton class to cache database objects.
Definition: DBStore.h:42
Belle2::StoreObjPtr::construct
bool construct(Args &&... params)
Construct an object of type T in this StoreObjPtr, using the provided constructor arguments.
Definition: StoreObjPtr.h:128
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::StoreObjPtr
Type-safe access to single objects in the data store.
Definition: ParticleList.h:33
Belle2::CalibrationAlgorithm::c_Failure
@ c_Failure
Failed =3 in Python.
Definition: CalibrationAlgorithm.h:54
Belle2::DBStore::Instance
static DBStore & Instance()
Instance of a singleton DBStore.
Definition: DBStore.cc:36
Belle2::ECLReferenceCrystalPerCrateCalib
General DB object to store one reference crystal per per ECL crate for calibration purposes.
Definition: ECLReferenceCrystalPerCrateCalib.h:51
Belle2::ECL::EclConfiguration::m_rf
static constexpr double m_rf
accelerating RF, http://ptep.oxfordjournals.org/content/2013/3/03A006.full.pdf
Definition: EclConfiguration.h:45
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::DBStore::updateEvent
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:150
Belle2::CalibrationAlgorithm::getRunList
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Definition: CalibrationAlgorithm.h:276
Belle2::ECL::eclBhabhaTAlgorithm::cellIDLo
int cellIDLo
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
Definition: eclBhabhaTAlgorithm.h:50
Belle2::ECL::eclBhabhaTAlgorithm::meanCleanRebinFactor
double meanCleanRebinFactor
Rebinning factor for mean calculation.
Definition: eclBhabhaTAlgorithm.h:52
Belle2::ECL::ECLChannelMapper::initFromDB
bool initFromDB()
Initialize channel mapper from the conditions database.
Definition: ECLChannelMapper.cc:105
Belle2::ECLCrystalCalib::setCalibVector
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
Definition: ECLCrystalCalib.h:48
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::ECL::eclBhabhaTAlgorithm::debugFilenameBase
std::string debugFilenameBase
Name of file with debug output, eclBhabhaTAlgorithm.root by default.
Definition: eclBhabhaTAlgorithm.h:62
Belle2::ECL::eclBhabhaTAlgorithm::meanCleanCutMinFactor
double meanCleanCutMinFactor
After rebinning, create a mask for bins that have values less than meanCleanCutMinFactor times the ma...
Definition: eclBhabhaTAlgorithm.h:53
Belle2::DBStore::update
void update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:87