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