Belle II Software  release-08-00-10
eclTValidationAlgorithm Class Reference

Validate the ecl timing calibrations using a hadronic event selection. More...

#include <eclTValidationAlgorithm.h>

Inheritance diagram for eclTValidationAlgorithm:
Collaboration diagram for eclTValidationAlgorithm:

Public Types

enum  EResult {
  c_OK ,
  c_Iterate ,
  c_NotEnoughData ,
  c_Failure ,
  c_Undefined
}
 The result of calibration. More...
 

Public Member Functions

 eclTValidationAlgorithm ()
 ..Constructor
 
 eclTValidationAlgorithm (std::string physicsProcessCollectorName)
 ..Constructor - main one as it allows user to choose which collector data to analyse
 
 ~eclTValidationAlgorithm ()
 ..Destructor
 
std::string getPrefix () const
 Get the prefix used for getting calibration data.
 
bool checkPyExpRun (PyObject *pyObj)
 Checks that a PyObject can be successfully converted to an ExpRun type. More...
 
Calibration::ExpRun convertPyExpRun (PyObject *pyObj)
 Performs the conversion of PyObject to ExpRun. More...
 
std::string getCollectorName () const
 Alias for prefix. More...
 
void setPrefix (const std::string &prefix)
 Set the prefix used to identify datastore objects.
 
void setInputFileNames (PyObject *inputFileNames)
 Set the input file names used for this algorithm from a Python list. More...
 
PyObject * getInputFileNames ()
 Get the input file names used for this algorithm and pass them out as a Python list of unicode strings.
 
std::vector< Calibration::ExpRun > getRunListFromAllData () const
 Get the complete list of runs from inspection of collected data.
 
RunRange getRunRangeFromAllData () const
 Get the complete RunRange from inspection of collected data.
 
IntervalOfValidity getIovFromAllData () const
 Get the complete IoV from inspection of collected data.
 
void fillRunToInputFilesMap ()
 Fill the mapping of ExpRun -> Files.
 
std::string getGranularity () const
 Get the granularity of collected data.
 
EResult execute (std::vector< Calibration::ExpRun > runs={}, int iteration=0, IntervalOfValidity iov=IntervalOfValidity())
 Runs calibration over vector of runs for a given iteration. More...
 
EResult execute (PyObject *runs, int iteration=0, IntervalOfValidity iov=IntervalOfValidity())
 Runs calibration over Python list of runs. Converts to C++ and then calls the other execute() function.
 
std::list< Database::DBImportQuery > & getPayloads ()
 Get constants (in TObjects) for database update from last execution.
 
std::list< Database::DBImportQuerygetPayloadValues ()
 Get constants (in TObjects) for database update from last execution but passed by VALUE.
 
bool commit ()
 Submit constants from last calibration into database.
 
bool commit (std::list< Database::DBImportQuery > payloads)
 Submit constants from a (potentially previous) set of payloads.
 
const std::string & getDescription () const
 Get the description of the algoithm (set by developers in constructor)
 
bool loadInputJson (const std::string &jsonString)
 Load the m_inputJson variable from a string (useful from Python interface). The rturn bool indicates success or failure.
 
const std::string dumpOutputJson () const
 Dump the JSON string of the output JSON object.
 
const std::vector< Calibration::ExpRun > findPayloadBoundaries (std::vector< Calibration::ExpRun > runs, int iteration=0)
 Used to discover the ExpRun boundaries that you want the Python CAF to execute on. This is optional and only used in some.
 
template<>
std::shared_ptr< TTree > getObjectPtr (const std::string &name, const std::vector< Calibration::ExpRun > &requestedRuns)
 Specialization of getObjectPtr<TTree>.
 

Public Attributes

int cellIDLo
 Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
 
int cellIDHi
 Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
 
bool readPrevCrysPayload
 Read the previous crystal payload values for comparison.
 
double meanCleanRebinFactor
 Rebinning factor for mean calculation.
 
double meanCleanCutMinFactor
 After rebinning, create a mask for bins that have values less than meanCleanCutMinFactor times the maximum bin value. More...
 
double clusterTimesFractionWindow_maxtime
 Maximum time for window to calculate cluster time fraction, in ns.
 
bool debugOutput
 Save every histogram and fitted function to debugFilename.
 
std::string debugFilenameBase
 Name of file with debug output, eclTValidationAlgorithm.root by default.
 

Protected Member Functions

EResult calibrate () override
 ..Run algorithm on events More...
 
void setInputFileNames (std::vector< std::string > inputFileNames)
 Set the input file names used for this algorithm. More...
 
virtual bool isBoundaryRequired (const Calibration::ExpRun &)
 Given the current collector data, make a decision about whether or not this run should be the start of a payload boundary.
 
virtual void boundaryFindingSetup (std::vector< Calibration::ExpRun >, int)
 If you need to make some changes to your algorithm class before 'findPayloadBoundaries' is run, make them in this function.
 
virtual void boundaryFindingTearDown ()
 Put your algorithm back into a state ready for normal execution if you need to.
 
const std::vector< Calibration::ExpRun > & getRunList () const
 Get the list of runs for which calibration is called.
 
int getIteration () const
 Get current iteration.
 
std::vector< std::string > getVecInputFileNames () const
 Get the input file names used for this algorithm as a STL vector.
 
template<class T >
std::shared_ptr< T > getObjectPtr (const std::string &name, const std::vector< Calibration::ExpRun > &requestedRuns)
 Get calibration data object by name and list of runs, the Merge function will be called to generate the overall object.
 
template<class T >
std::shared_ptr< T > getObjectPtr (std::string name)
 Get calibration data object (for all runs the calibration is requested for) This function will only work during or after execute() has been called once.
 
template<>
shared_ptr< TTree > getObjectPtr (const string &name, const vector< ExpRun > &requestedRuns)
 We cheekily cast the TChain to TTree for the returned pointer so that the user never knows Hopefully this doesn't cause issues if people do low level stuff to the tree...
 
std::string getGranularityFromData () const
 Get the granularity of collected data.
 
void saveCalibration (TClonesArray *data, const std::string &name)
 Store DBArray payload with given name with default IOV.
 
void saveCalibration (TClonesArray *data, const std::string &name, const IntervalOfValidity &iov)
 Store DBArray with given name and custom IOV.
 
void saveCalibration (TObject *data)
 Store DB payload with default name and default IOV.
 
void saveCalibration (TObject *data, const IntervalOfValidity &iov)
 Store DB payload with default name and custom IOV.
 
void saveCalibration (TObject *data, const std::string &name)
 Store DB payload with given name with default IOV.
 
void saveCalibration (TObject *data, const std::string &name, const IntervalOfValidity &iov)
 Store DB payload with given name and custom 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)
 
void clearCalibrationData ()
 Clear calibration data.
 
Calibration::ExpRun getAllGranularityExpRun () const
 Returns the Exp,Run pair that means 'Everything'. Currently unused.
 
void resetInputJson ()
 Clears the m_inputJson member variable.
 
void resetOutputJson ()
 Clears the m_outputJson member variable.
 
template<class T >
void setOutputJsonValue (const std::string &key, const T &value)
 Set a key:value pair for the outputJson object, expected to used interally during calibrate()
 
template<class T >
const T getOutputJsonValue (const std::string &key) const
 Get a value using a key from the JSON output object, not sure why you would want to do this.
 
template<class T >
const T getInputJsonValue (const std::string &key) const
 Get an input JSON value using a key. The normal exceptions are raised when the key doesn't exist.
 
const nlohmann::json & getInputJsonObject () const
 Get the entire top level JSON object. We explicitly say this must be of object type so that we might pick.
 
bool inputJsonKeyExists (const std::string &key) const
 Test for a key in the input JSON object.
 

Protected Attributes

std::vector< Calibration::ExpRun > m_boundaries
 When using the boundaries functionality from isBoundaryRequired, this is used to store the boundaries. It is cleared when.
 

Private Member Functions

std::string getExpRunString (Calibration::ExpRun &expRun) const
 Gets the "exp.run" string repr. of (exp,run)
 
std::string getFullObjectPath (const std::string &name, Calibration::ExpRun expRun) const
 constructs the full TDirectory + Key name of an object in a TFile based on its name and exprun
 

Private Attributes

std::vector< std::string > m_inputFileNames
 List of input files to the Algorithm, will initially be user defined but then gets the wildcards expanded during execute()
 
std::map< Calibration::ExpRun, std::vector< std::string > > m_runsToInputFiles
 Map of Runs to input files. Gets filled when you call getRunRangeFromAllData, gets cleared when setting input files again.
 
std::string m_granularityOfData
 Granularity of input data. This only changes when the input files change so it isn't specific to an execution.
 
ExecutionData m_data
 Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.
 
std::string m_description {""}
 Description of the algorithm.
 
std::string m_prefix {""}
 The name of the TDirectory the collector objects are contained within.
 
nlohmann::json m_jsonExecutionInput = nlohmann::json::object()
 Optional input JSON object used to make decisions about how to execute the algorithm code.
 
nlohmann::json m_jsonExecutionOutput = nlohmann::json::object()
 Optional output JSON object that can be set during the execution by the underlying algorithm code.
 

Static Private Attributes

static const Calibration::ExpRun m_allExpRun = make_pair(-1, -1)
 allExpRun
 

Detailed Description

Validate the ecl timing calibrations using a hadronic event selection.

Definition at line 31 of file eclTValidationAlgorithm.h.

Member Enumeration Documentation

◆ EResult

enum EResult
inherited

The result of calibration.

Enumerator
c_OK 

Finished successfuly =0 in Python.

c_Iterate 

Needs iteration =1 in Python.

c_NotEnoughData 

Needs more data =2 in Python.

c_Failure 

Failed =3 in Python.

c_Undefined 

Not yet known (before execution) =4 in Python.

Definition at line 40 of file CalibrationAlgorithm.h.

Member Function Documentation

◆ calibrate()

CalibrationAlgorithm::EResult calibrate ( )
overrideprotectedvirtual

..Run algorithm on events

Put root into batch mode so that we don't try to open a graphics window

Write out job parameters


Implements CalibrationAlgorithm.

Definition at line 87 of file eclTValidationAlgorithm.cc.

88 {
90  gROOT->SetBatch();
91 
92 
94  B2INFO("eclTValidationAlgorithm parameters:");
95  B2INFO("cellIDLo = " << cellIDLo);
96  B2INFO("cellIDHi = " << cellIDHi);
97  B2INFO("readPrevCrysPayload = " << readPrevCrysPayload);
98  B2INFO("meanCleanRebinFactor = " << meanCleanRebinFactor);
99  B2INFO("meanCleanCutMinFactor = " << meanCleanCutMinFactor);
100  B2INFO("clusterTimesFractionWindow_maxtime = " << clusterTimesFractionWindow_maxtime);
101 
102 
103  /* Histogram with the data collected by eclTimeCalibrationValidationCollector*/
104  auto clusterTime = getObjectPtr<TH1F>("clusterTime");
105  auto clusterTime_cid = getObjectPtr<TH2F>("clusterTime_cid");
106  auto clusterTime_run = getObjectPtr<TH2F>("clusterTime_run");
107  auto clusterTimeClusterE = getObjectPtr<TH2F>("clusterTimeClusterE");
108  auto dt99_clusterE = getObjectPtr<TH2F>("dt99_clusterE");
109  auto eventT0 = getObjectPtr<TH1F>("eventT0");
110  auto clusterTimeE0E1diff = getObjectPtr<TH1F>("clusterTimeE0E1diff");
111 
112  // Collect other plots just for reference - combines all the runs for these plots.
113  auto cutflow = getObjectPtr<TH1F>("cutflow");
114 
115  vector <int> binProjectionLeft_Time_vs_E_runDep ;
116  vector <int> binProjectionRight_Time_vs_E_runDep ;
117 
118  for (int binCounter = 1; binCounter <= clusterTimeClusterE->GetNbinsX(); binCounter++) {
119  binProjectionLeft_Time_vs_E_runDep.push_back(binCounter);
120  binProjectionRight_Time_vs_E_runDep.push_back(binCounter);
121  }
122 
123  if (!clusterTime_cid) return c_Failure;
124 
127  TFile* histfile = 0;
128 
129  // Vector of time offsets to track how far from nominal the cluster times are.
130  vector<float> t_offsets(ECLElementNumbers::c_NCrystals, 0.0);
131  vector<float> t_offsets_unc(ECLElementNumbers::c_NCrystals, 0.0);
132  vector<long> numClusterPerCrys(ECLElementNumbers::c_NCrystals, 0);
133  vector<bool> crysHasGoodFitandStats(ECLElementNumbers::c_NCrystals, false);
134  vector<bool> crysHasGoodFit(ECLElementNumbers::c_NCrystals, false);
135  int numCrysWithNonZeroEntries = 0 ;
136  int numCrysWithGoodFit = 0 ;
137 
138  int minNumEntries = 40;
139 
140  double mean;
141  double sigma;
142 
143 
144  bool minRunNumBool = false;
145  bool maxRunNumBool = false;
146  int minRunNum = -1;
147  int maxRunNum = -1;
148  int minExpNum = -1;
149  int maxExpNum = -1;
150  for (auto expRun : getRunList()) {
151  int expNumber = expRun.first;
152  int runNumber = expRun.second;
153  if (!minRunNumBool) {
154  minExpNum = expNumber;
155  minRunNum = runNumber;
156  minRunNumBool = true;
157  }
158  if (!maxRunNumBool) {
159  maxExpNum = expNumber;
160  maxRunNum = runNumber;
161  maxRunNumBool = true;
162  }
163  if (((minRunNum > runNumber) && (minExpNum >= expNumber)) ||
164  (minExpNum > expNumber)) {
165  minExpNum = expNumber;
166  minRunNum = runNumber;
167  }
168  if (((maxRunNum < runNumber) && (maxExpNum <= expNumber)) ||
169  (maxExpNum < expNumber))
170 
171  {
172  maxExpNum = expNumber;
173  maxRunNum = runNumber;
174  }
175  }
176 
177  B2INFO("debugFilenameBase = " << debugFilenameBase);
178  string runNumsString = string("_") + to_string(minExpNum) + "_" + to_string(minRunNum) + string("-") +
179  to_string(maxExpNum) + "_" + to_string(maxRunNum);
180  string debugFilename = debugFilenameBase + runNumsString + string(".root");
181 
182 
183  // Need to load information about the event/run/experiment to get the right database information
184  // Will be used for:
185  // * ECLChannelMapper (to map crystal to crates)
186  // * crystal payload updating for iterating crystal and crate fits
187  int eventNumberForCrates = 1;
188 
189 
190  //-------------------------------------------------------------------
191  /* Uploading older payloads for the current set of runs */
192 
194  // simulate the initialize() phase where we can register objects in the DataStore
196  evtPtr.registerInDataStore();
198  // now construct the event metadata
199  evtPtr.construct(eventNumberForCrates, minRunNum, minExpNum);
200  // and update the database contents
201  DBStore& dbstore = DBStore::Instance();
202  dbstore.update();
203  // this is only needed it the payload might be intra-run dependent,
204  // that is if it might change during one run as well
205  dbstore.updateEvent();
206 
207 
208  B2INFO("Uploading payload for exp " << minExpNum << ", run " << minRunNum << ", event " << eventNumberForCrates);
209  updateDBObjPtrs(eventNumberForCrates, minRunNum, minExpNum);
210  unique_ptr<ECLChannelMapper> crystalMapper(new ECL::ECLChannelMapper());
211  crystalMapper->initFromDB();
212 
213  /* 1/(4fRF) = 0.4913 ns/clock tick, where fRF is the accelerator RF frequency.
214  Same for all crystals. */
215  const double TICKS_TO_NS = 1.0 / (4.0 * EclConfiguration::getRF()) * 1e3;
216 
217  //------------------------------------------------------------------------
218  //..Read payloads from database
219  DBObjPtr<Belle2::ECLCrystalCalib> crystalTimeObject("ECLCrystalTimeOffset");
220  B2INFO("Dumping payload");
221 
222  //..Get vectors of values from the payloads
223  std::vector<float> currentValuesCrys = crystalTimeObject->getCalibVector();
224  std::vector<float> currentUncCrys = crystalTimeObject->getCalibUncVector();
225 
226  //..Print out a few values for quality control
227  B2INFO("Values read from database. Write out for their values for comparison against those from tcol");
228  for (int ic = 0; ic < ECLElementNumbers::c_NCrystals; ic += 500) {
229  B2INFO("ts: cellID " << ic + 1 << " " << currentValuesCrys[ic] << " +/- " << currentUncCrys[ic]);
230  }
231 
232 
233  //..Read in the previous crystal payload values
234  DBObjPtr<Belle2::ECLCrystalCalib> customPrevCrystalTimeObject("ECLCrystalTimeOffsetPreviousValues");
235  vector<float> prevValuesCrys(ECLElementNumbers::c_NCrystals);
236  if (readPrevCrysPayload) {
237  //..Get vectors of values from the payloads
238  prevValuesCrys = customPrevCrystalTimeObject->getCalibVector();
239 
240  //..Print out a few values for quality control
241  B2INFO("Previous values read from database. Write out for their values for comparison against those from tcol");
242  for (int ic = 0; ic < ECLElementNumbers::c_NCrystals; ic += 500) {
243  B2INFO("ts custom previous payload: cellID " << ic + 1 << " " << prevValuesCrys[ic]);
244  }
245  }
246 
247 
248  //------------------------------------------------------------------------
249  //..Start looking at timing information
250 
251  B2INFO("Debug output rootfile: " << debugFilename);
252  histfile = new TFile(debugFilename.c_str(), "recreate");
253 
254 
255  clusterTime ->Write();
256  clusterTime_cid ->Write();
257  clusterTime_run ->Write();
258  clusterTimeClusterE ->Write();
259  dt99_clusterE ->Write();
260  eventT0 ->Write();
261  clusterTimeE0E1diff ->Write();
262 
263  cutflow->Write();
264 
265 
266  double hist_tmin = clusterTime->GetXaxis()->GetXmin();
267  double hist_tmax = clusterTime->GetXaxis()->GetXmax();
268  int hist_nTbins = clusterTime->GetNbinsX();
269 
270  B2INFO("hist_tmin = " << hist_tmin);
271  B2INFO("hist_tmax = " << hist_tmax);
272  B2INFO("hist_nTbins = " << hist_nTbins);
273 
274  double time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
275  double time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
276 
277 
278  // define histogram for keeping track of the peak of the cluster times per crystal
279  auto peakClusterTime_cid = new TH1F("peakClusterTime_cid", ";cell id;Peak cluster time [ns]", ECLElementNumbers::c_NCrystals, 1,
281  auto peakClusterTimes = new TH1F("peakClusterTimes",
282  "-For crystals with at least one hit-;Peak cluster time [ns];Number of crystals",
283  hist_nTbins, hist_tmin, hist_tmax);
284  auto peakClusterTimesGoodFit = new TH1F("peakClusterTimesGoodFit",
285  "-For crystals with a good fit to distribution of hits-;Peak cluster time [ns];Number of crystals",
286  hist_nTbins, hist_tmin, hist_tmax);
287 
288  auto peakClusterTimesGoodFit__cid = new TH1F("peakClusterTimesGoodFit__cid",
289  "-For crystals with a good fit to distribution of hits-;cell id (only crystals with good fit);Peak cluster time [ns]",
291 
292 
293  // define histograms to keep track of the difference in the new crystal times vs the old ones
294  auto tsNew_MINUS_tsCustomPrev__cid = new TH1F("TsNew_MINUS_TsCustomPrev__cid",
295  ";cell id; ts(new|merged) - ts(old = 'pre-calib'|merged) [ns]",
297 
298  auto tsNew_MINUS_tsCustomPrev = new TH1F("TsNew_MINUS_TsCustomPrev",
299  ";ts(new | merged) - ts(old = 'pre-calib' | merged) [ns];Number of crystals",
300  285, -69.5801, 69.5801);
301 
302 
303 
304  // Histogram to keep track of the fraction of cluster times within a window.
305  double timeWindow_maxTime = clusterTimesFractionWindow_maxtime;
306  B2INFO("timeWindow_maxTime = " << timeWindow_maxTime);
307  int binyLeft = clusterTime_cid->GetYaxis()->FindBin(-timeWindow_maxTime);
308  int binyRight = clusterTime_cid->GetYaxis()->FindBin(timeWindow_maxTime);
309  double windowLowTimeFromBin = clusterTime_cid->GetYaxis()->GetBinLowEdge(binyLeft);
310  double windowHighTimeFromBin = clusterTime_cid->GetYaxis()->GetBinLowEdge(binyRight + 1);
311  std::string s_lowTime = std::to_string(windowLowTimeFromBin);
312  std::string s_highTime = std::to_string(windowHighTimeFromBin);
313  TString fracWindowTitle = "Fraction of cluster times in window [" + s_lowTime + ", " + s_highTime +
314  "] ns;cell id;Fraction of cluster times in window";
315  B2INFO("fracWindowTitle = " << fracWindowTitle);
316  TString fracWindowInGoodECLRingsTitle = "Fraction of cluster times in window [" + s_lowTime + ", " + s_highTime +
317  "] ns and in good ECL theta rings;cell id;Fraction cluster times in window + good ECL rings";
318  B2INFO("fracWindowInGoodECLRingsTitle = " << fracWindowInGoodECLRingsTitle);
319  B2INFO("Good ECL rings skip gaps in the acceptance, and includes ECL theta IDs: 3-10, 15-39, 44-56, 61-66.");
320 
321  TString fracWindowHistTitle = "Fraction of cluster times in window [" + s_lowTime + ", " + s_highTime +
322  "] ns;Fraction of cluster times in window;Number of crystals";
323 
324  auto clusterTimeNumberInWindow__cid = new TH1F("clusterTimeNumberInWindow__cid", fracWindowTitle, ECLElementNumbers::c_NCrystals, 1,
326  auto clusterTimeNumberInWindowInGoodECLRings__cid = new TH1F("clusterTimeNumberInWindowInGoodECLRings__cid", fracWindowTitle,
329  auto clusterTimeNumber__cid = new TH1F("clusterTimeNumber_cid", fracWindowTitle, ECLElementNumbers::c_NCrystals, 1,
331  auto clusterTimeFractionInWindow = new TH1F("clusterTimeFractionInWindow", fracWindowHistTitle, 110, 0.0, 1.1);
332 
333  clusterTimeNumberInWindow__cid->Sumw2();
334  clusterTimeNumberInWindowInGoodECLRings__cid->Sumw2();
335  clusterTimeNumber__cid->Sumw2();
336 
337 
338 
339  /* CRYSTAL BY CRYSTAL VALIDATION */
340 
342 
343  // Loop over all the crystals for doing the crystal calibation
344  for (int crys_id = cellIDLo; crys_id <= cellIDHi; crys_id++) {
345  double clusterTime_mean = 0;
346  double clusterTime_mean_unc = 0;
347 
348  B2INFO("Crystal cell id = " << crys_id);
349 
350  eclgeo->Mapping(crys_id - 1);
351  int thetaID = eclgeo->GetThetaID();
352 
353 
354  /* Determining which bins to mask out for mean calculation
355  */
356 
357  TH1D* h_time = clusterTime_cid->ProjectionY((std::string("h_time_psi__") + std::to_string(crys_id)).c_str(),
358  crys_id, crys_id);
359  TH1D* h_timeMask = (TH1D*)h_time->Clone();
360  TH1D* h_timeMasked = (TH1D*)h_time->Clone((std::string("h_time_psi_masked__") + std::to_string(crys_id)).c_str());
361  TH1D* h_timeRebin = (TH1D*)h_time->Clone();
362 
363  // Do rebinning and cleaning of some bins but only if the user selection values call for it since it slows the code down
364  if (meanCleanRebinFactor != 1 || meanCleanCutMinFactor != 1) {
365 
366  h_timeRebin->Rebin(meanCleanRebinFactor);
367 
368  h_timeMask->Scale(0.0); // set all bins to being masked off
369 
370  time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
371  time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
372 
373  // Find value of bin with max value
374  double histRebin_max = h_timeRebin->GetMaximum();
375 
376  bool maskedOutNonZeroBin = false;
377  // Loop over all bins to find those with content less than a certain threshold. Mask the non-rebinned histogram for the corresponding bins
378  for (int bin = 1; bin <= h_timeRebin->GetNbinsX(); bin++) {
379  for (int rebinCounter = 1; rebinCounter <= meanCleanRebinFactor; rebinCounter++) {
380  int nonRebinnedBinNumber = (bin - 1) * meanCleanRebinFactor + rebinCounter;
381  if (nonRebinnedBinNumber < h_time->GetNbinsX()) {
382  if (h_timeRebin->GetBinContent(bin) >= histRebin_max * meanCleanCutMinFactor) {
383  h_timeMask->SetBinContent(nonRebinnedBinNumber, 1);
384 
385  // Save the lower and upper edges of the rebin histogram time range for fitting purposes
386  double x_lower = h_timeRebin->GetXaxis()->GetBinLowEdge(bin);
387  double x_upper = h_timeRebin->GetXaxis()->GetBinUpEdge(bin);
388  if (x_lower < time_fit_min) {
389  time_fit_min = x_lower;
390  }
391  if (x_upper > time_fit_max) {
392  time_fit_max = x_upper;
393  }
394 
395  } else {
396  if (h_time->GetBinContent(nonRebinnedBinNumber) > 0) {
397  B2DEBUG(22, "Setting bin " << nonRebinnedBinNumber << " from " << h_timeMasked->GetBinContent(nonRebinnedBinNumber) << " to 0");
398  maskedOutNonZeroBin = true;
399  }
400  h_timeMasked->SetBinContent(nonRebinnedBinNumber, 0);
401  }
402  }
403  }
404  }
405  B2INFO("Bins with non-zero values have been masked out: " << maskedOutNonZeroBin);
406  h_timeMasked->ResetStats();
407  h_timeMask->ResetStats();
408 
409  }
410 
411  // Calculate mean from masked histogram
412  double default_meanMasked = h_timeMasked->GetMean();
413  //double default_meanMasked_unc = h_timeMasked->GetMeanError();
414  B2INFO("default_meanMasked = " << default_meanMasked);
415 
416 
417  // Get the overall mean and standard deviation of the distribution within the plot. This doesn't require a fit.
418  double default_mean = h_time->GetMean();
419  double default_mean_unc = h_time->GetMeanError();
420  double default_sigma = h_time->GetStdDev();
421 
422  B2INFO("Fitting crystal between " << time_fit_min << " and " << time_fit_max);
423 
424  // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
425  TF1* gaus = new TF1("func", "gaus(0)", time_fit_min, time_fit_max);
426  gaus->SetParNames("numCrystalHitsNormalization", "mean", "sigma");
427  /*
428  gaus->ReleaseParameter(0); // number of crystals
429  gaus->ReleaseParameter(1); // mean
430  gaus->ReleaseParameter(2); // standard deviation
431  */
432 
433  double hist_max = h_time->GetMaximum();
434 
435  //=== Estimate initial value of sigma as std dev.
436  double stddev = h_time->GetStdDev();
437  sigma = stddev;
438  mean = default_mean;
439 
440  //=== Setting parameters for initial iteration
441  gaus->SetParameter(0, hist_max / 2.);
442  gaus->SetParameter(1, mean);
443  gaus->SetParameter(2, sigma);
444  // L -- Use log likelihood method
445  // I -- Use integral of function in bin instead of value at bin center // not using
446  // R -- Use the range specified in the function range
447  // B -- Fix one or more parameters with predefined function // not using
448  // Q -- Quiet mode
449 
450  h_timeMasked->Fit(gaus, "LQR"); // L for likelihood, R for x-range, Q for fit quiet mode
451 
452  double fit_mean = gaus->GetParameter(1);
453  double fit_mean_unc = gaus->GetParError(1);
454  double fit_sigma = gaus->GetParameter(2);
455 
456  double meanDiff = fit_mean - default_mean;
457  double meanUncDiff = fit_mean_unc - default_mean_unc;
458  double sigmaDiff = fit_sigma - default_sigma;
459 
460  bool good_fit = false;
461 
462  if ((fabs(meanDiff) > 10) ||
463  (fabs(meanUncDiff) > 10) ||
464  (fabs(sigmaDiff) > 10) ||
465  (fit_mean_unc > 3) ||
466  (fit_sigma < 0.1) ||
467  (fit_mean < time_fit_min) ||
468  (fit_mean > time_fit_max)) {
469  B2INFO("Crystal cell id = " << crys_id);
470  B2INFO("fit mean, default mean = " << fit_mean << ", " << default_mean);
471  B2INFO("fit mean unc, default mean unc = " << fit_mean_unc << ", " << default_mean_unc);
472  B2INFO("fit sigma, default sigma = " << fit_sigma << ", " << default_sigma);
473 
474  B2INFO("crystal fit mean - hist mean = " << meanDiff);
475  B2INFO("fit mean unc. - hist mean unc. = " << meanUncDiff);
476  B2INFO("fit sigma - hist sigma = " << sigmaDiff);
477 
478  B2INFO("fit_mean = " << fit_mean);
479  B2INFO("time_fit_min = " << time_fit_min);
480  B2INFO("time_fit_max = " << time_fit_max);
481 
482  if (fabs(meanDiff) > 10) B2INFO("fit mean diff too large");
483  if (fabs(meanUncDiff) > 10) B2INFO("fit mean unc diff too large");
484  if (fabs(sigmaDiff) > 10) B2INFO("fit mean sigma diff too large");
485  if (fit_mean_unc > 3) B2INFO("fit mean unc too large");
486  if (fit_sigma < 0.1) B2INFO("fit sigma too small");
487 
488  } else {
489  good_fit = true;
490  numCrysWithGoodFit++;
491  crysHasGoodFit[crys_id - 1] = true ;
492  }
493 
494 
495  int numEntries = h_time->GetEntries();
496  // If number of entries in histogram is greater than X then use the statistical information from the data otherwise leave crystal uncalibrated. Histograms are still shown though.
497  // ALSO require the that fits are good.
498  if ((numEntries >= minNumEntries) && good_fit) {
499  clusterTime_mean = fit_mean;
500  clusterTime_mean_unc = fit_mean_unc;
501  crysHasGoodFitandStats[crys_id - 1] = true ;
502  } else {
503  clusterTime_mean = default_mean;
504  clusterTime_mean_unc = default_mean_unc;
505  }
506 
507  if (numEntries < minNumEntries) B2INFO("Number of entries less than minimum");
508  if (numEntries == 0) B2INFO("Number of entries == 0");
509 
510 
511  t_offsets[crys_id - 1] = clusterTime_mean ;
512  t_offsets_unc[crys_id - 1] = clusterTime_mean_unc ;
513  numClusterPerCrys[crys_id - 1] = numEntries ;
514 
515  histfile->WriteTObject(h_time, (std::string("h_time_psi") + std::to_string(crys_id)).c_str());
516  histfile->WriteTObject(h_timeMasked, (std::string("h_time_psi_masked") + std::to_string(crys_id)).c_str());
517 
518  // Set this for each crystal even if there are zero entries
519  peakClusterTime_cid->SetBinContent(crys_id, t_offsets[crys_id - 1]);
520  peakClusterTime_cid->SetBinError(crys_id, t_offsets_unc[crys_id - 1]);
521 
522  /* Store mean cluster time info in a separate histogram but only if there is at
523  least one entry for that crystal. */
524  if (numEntries > 0) {
525  peakClusterTimes->Fill(t_offsets[crys_id - 1]);
526  numCrysWithNonZeroEntries++ ;
527  }
528  if ((numEntries >= minNumEntries) && good_fit) {
529  peakClusterTimesGoodFit->Fill(t_offsets[crys_id - 1]);
530  peakClusterTimesGoodFit__cid->SetBinContent(crys_id, t_offsets[crys_id - 1]);
531  peakClusterTimesGoodFit__cid->SetBinError(crys_id, t_offsets_unc[crys_id - 1]);
532  }
533 
534 
535  // Find the fraction of cluster times within +-X ns and fill histograms
536  double numClusterTimesWithinWindowFraction = h_time->Integral(binyLeft, binyRight) ;
537  double clusterTimesWithinWindowFraction = numClusterTimesWithinWindowFraction;
538  if (numEntries > 0) {
539  clusterTimesWithinWindowFraction /= numEntries;
540  } else {
541  clusterTimesWithinWindowFraction = -0.1;
542  }
543 
544  B2INFO("Crystal cell id = " << crys_id << ", theta id = " <<
545  thetaID << ", clusterTimesWithinWindowFraction = " <<
546  numClusterTimesWithinWindowFraction << " / " << numEntries << " = " <<
547  clusterTimesWithinWindowFraction);
548 
549  clusterTimeFractionInWindow->Fill(clusterTimesWithinWindowFraction);
550  clusterTimeNumberInWindow__cid->SetBinContent(crys_id, numClusterTimesWithinWindowFraction);
551  clusterTimeNumber__cid->SetBinContent(crys_id, numEntries);
552 
553  if ((thetaID >= 3 && thetaID <= 10) ||
554  (thetaID >= 15 && thetaID <= 39) ||
555  (thetaID >= 44 && thetaID <= 56) ||
556  (thetaID >= 61 && thetaID <= 66)) {
557  clusterTimeNumberInWindowInGoodECLRings__cid->SetBinContent(crys_id, numClusterTimesWithinWindowFraction);
558  }
559 
560 
561  delete gaus;
562  }
563 
564  // Find the fraction of cluster times within +-X ns and fill histogram
565  auto g_clusterTimeFractionInWindow__cid = new TGraphAsymmErrors(clusterTimeNumberInWindow__cid, clusterTimeNumber__cid, "w");
566  auto g_clusterTimeFractionInWindowInGoodECLRings__cid = new TGraphAsymmErrors(clusterTimeNumberInWindowInGoodECLRings__cid,
567  clusterTimeNumber__cid, "w");
568  g_clusterTimeFractionInWindow__cid->SetTitle(fracWindowTitle);
569  g_clusterTimeFractionInWindowInGoodECLRings__cid->SetTitle(fracWindowInGoodECLRingsTitle);
570 
571 
572  peakClusterTime_cid->ResetStats();
573  peakClusterTimesGoodFit__cid->ResetStats();
574 
575  histfile->WriteTObject(peakClusterTime_cid, "peakClusterTime_cid");
576  histfile->WriteTObject(peakClusterTimes, "peakClusterTimes");
577  histfile->WriteTObject(peakClusterTimesGoodFit__cid, "peakClusterTimesGoodFit__cid");
578  histfile->WriteTObject(peakClusterTimesGoodFit, "peakClusterTimesGoodFit");
579  histfile->WriteTObject(g_clusterTimeFractionInWindow__cid, "g_clusterTimeFractionInWindow__cid");
580  histfile->WriteTObject(g_clusterTimeFractionInWindowInGoodECLRings__cid, "g_clusterTimeFractionInWindowInGoodECLRings__cid");
581  histfile->WriteTObject(clusterTimeFractionInWindow, "clusterTimeFractionInWindow");
582 
583 
584 
585  /* -----------------------------------------------------------
586  Fit the time histograms for different energy slices */
587 
588  vector <int> binProjectionLeft = binProjectionLeft_Time_vs_E_runDep;
589  vector <int> binProjectionRight = binProjectionRight_Time_vs_E_runDep;
590 
591  auto h2 = clusterTimeClusterE;
592 
593 
594  double max_E = h2->GetXaxis()->GetXmax();
595 
596  // Determine the energy bins. Save the left edge for histogram purposes
597  vector <double> E_binEdges(binProjectionLeft.size() + 1);
598  for (long unsigned int x_bin = 0; x_bin < binProjectionLeft.size(); x_bin++) {
599  TH1D* h_E_t_slice = h2->ProjectionX("h_E_t_slice", 1, 1) ;
600  E_binEdges[x_bin] = h_E_t_slice->GetXaxis()->GetBinLowEdge(binProjectionLeft[x_bin]) ;
601  B2INFO("E_binEdges[" << x_bin << "] = " << E_binEdges[x_bin]);
602  if (x_bin == binProjectionLeft.size() - 1) {
603  E_binEdges[x_bin + 1] = max_E ;
604  B2INFO("E_binEdges[" << x_bin + 1 << "] = " << E_binEdges[x_bin + 1]);
605  }
606  }
607 
608 
609  auto clusterTimePeak_ClusterEnergy_varBin = new TH1F("clusterTimePeak_ClusterEnergy_varBin",
610  ";ECL cluster energy [GeV];Cluster time fit position [ns]", E_binEdges.size() - 1, &(E_binEdges[0]));
611  auto clusterTimePeakWidth_ClusterEnergy_varBin = new TH1F("clusterTimePeakWidth_ClusterEnergy_varBin",
612  ";ECL cluster energy [GeV];Cluster time fit width [ns]", E_binEdges.size() - 1, &(E_binEdges[0]));
613 
614  int Ebin_counter = 1 ;
615 
616  // Loop over all the energy bins
617  for (long unsigned int x_bin = 0; x_bin < binProjectionLeft.size(); x_bin++) {
618  double clusterTime_mean = 0;
619  double clusterTime_mean_unc = 0;
620  double clusterTime_sigma = 0;
621 
622  B2INFO("x_bin = " << x_bin);
623 
624  /* Determining which bins to mask out for mean calculation
625  */
626  TH1D* h_time = h2->ProjectionY(("h_time_E_slice_" + std::to_string(x_bin)).c_str(), binProjectionLeft[x_bin],
627  binProjectionRight[x_bin]) ;
628 
629 
630  TH1D* h_E_t_slice = h2->ProjectionX("h_E_t_slice", 1, 1) ;
631  double lowE = h_E_t_slice->GetXaxis()->GetBinLowEdge(binProjectionLeft[x_bin]) ;
632  double highE = h_E_t_slice->GetXaxis()->GetBinUpEdge(binProjectionRight[x_bin]) ;
633  double meanE = (lowE + highE) / 2.0 ;
634 
635  B2INFO("bin " << Ebin_counter << ": low E = " << lowE << ", high E = " << highE << " GeV");
636 
637  TH1D* h_timeMask = (TH1D*)h_time->Clone();
638  TH1D* h_timeMasked = (TH1D*)h_time->Clone((std::string("h_time_E_slice_masked__") + std::to_string(meanE)).c_str());
639  TH1D* h_timeRebin = (TH1D*)h_time->Clone();
640 
641 
642  if (meanCleanRebinFactor != 1 || meanCleanCutMinFactor != 1) {
643 
644  h_timeRebin->Rebin(meanCleanRebinFactor);
645 
646  h_timeMask->Scale(0.0); // set all bins to being masked off
647 
648  time_fit_min = hist_tmax; // Set min value to largest possible value so that it gets reset
649  time_fit_max = hist_tmin; // Set max value to smallest possible value so that it gets reset
650 
651  // Find value of bin with max value
652  double histRebin_max = h_timeRebin->GetMaximum();
653 
654  bool maskedOutNonZeroBin = false;
655  // Loop over all bins to find those with content less than a certain threshold. Mask the non-rebinned histogram for the corresponding bins
656  for (int bin = 1; bin <= h_timeRebin->GetNbinsX(); bin++) {
657  for (int rebinCounter = 1; rebinCounter <= meanCleanRebinFactor; rebinCounter++) {
658  int nonRebinnedBinNumber = (bin - 1) * meanCleanRebinFactor + rebinCounter;
659  if (nonRebinnedBinNumber < h_time->GetNbinsX()) {
660  if (h_timeRebin->GetBinContent(bin) >= histRebin_max * meanCleanCutMinFactor) {
661  h_timeMask->SetBinContent(nonRebinnedBinNumber, 1);
662 
663  // Save the lower and upper edges of the rebin histogram time range for fitting purposes
664  double x_lower = h_timeRebin->GetXaxis()->GetBinLowEdge(bin);
665  double x_upper = h_timeRebin->GetXaxis()->GetBinUpEdge(bin);
666  if (x_lower < time_fit_min) {
667  time_fit_min = x_lower;
668  }
669  if (x_upper > time_fit_max) {
670  time_fit_max = x_upper;
671  }
672 
673  } else {
674  if (h_time->GetBinContent(nonRebinnedBinNumber) > 0) {
675  B2DEBUG(22, "Setting bin " << nonRebinnedBinNumber << " from " << h_timeMasked->GetBinContent(nonRebinnedBinNumber) << " to 0");
676  maskedOutNonZeroBin = true;
677  }
678  h_timeMasked->SetBinContent(nonRebinnedBinNumber, 0);
679  }
680  }
681  }
682  }
683  B2INFO("Bins with non-zero values have been masked out: " << maskedOutNonZeroBin);
684  h_timeMasked->ResetStats();
685  h_timeMask->ResetStats();
686 
687  }
688 
689 
690  // Calculate mean from masked histogram
691  double default_meanMasked = h_timeMasked->GetMean();
692  //double default_meanMasked_unc = h_timeMasked->GetMeanError();
693  B2INFO("default_meanMasked = " << default_meanMasked);
694 
695 
696  // Get the overall mean and standard deviation of the distribution within the plot. This doesn't require a fit.
697  double default_mean = h_time->GetMean();
698  double default_mean_unc = h_time->GetMeanError();
699  double default_sigma = h_time->GetStdDev();
700 
701  B2INFO("Fitting crystal between " << time_fit_min << " and " << time_fit_max);
702 
703  // gaus(0) is a substitute for [0]*exp(-0.5*((x-[1])/[2])**2)
704  TF1* gaus = new TF1("func", "gaus(0)", time_fit_min, time_fit_max);
705  gaus->SetParNames("numCrystalHitsNormalization", "mean", "sigma");
706  /*
707  gaus->ReleaseParameter(0); // number of crystals
708  gaus->ReleaseParameter(1); // mean
709  gaus->ReleaseParameter(2); // standard deviation
710  */
711 
712  double hist_max = h_time->GetMaximum();
713 
714  //=== Estimate initial value of sigma as std dev.
715  double stddev = h_time->GetStdDev();
716  sigma = stddev;
717  mean = default_mean;
718 
719  //=== Setting parameters for initial iteration
720  gaus->SetParameter(0, hist_max / 2.);
721  gaus->SetParameter(1, mean);
722  gaus->SetParameter(2, sigma);
723  // L -- Use log likelihood method
724  // I -- Use integral of function in bin instead of value at bin center // not using
725  // R -- Use the range specified in the function range
726  // B -- Fix one or more parameters with predefined function // not using
727  // Q -- Quiet mode
728 
729  h_timeMasked->Fit(gaus, "LQR"); // L for likelihood, R for x-range, Q for fit quiet mode
730 
731  double fit_mean = gaus->GetParameter(1);
732  double fit_mean_unc = gaus->GetParError(1);
733  double fit_sigma = gaus->GetParameter(2);
734 
735  double meanDiff = fit_mean - default_mean;
736  double meanUncDiff = fit_mean_unc - default_mean_unc;
737  double sigmaDiff = fit_sigma - default_sigma;
738 
739  bool good_fit = false;
740 
741  if ((fabs(meanDiff) > 10) ||
742  (fabs(meanUncDiff) > 10) ||
743  (fabs(sigmaDiff) > 10) ||
744  (fit_mean_unc > 3) ||
745  (fit_sigma < 0.1) ||
746  (fit_mean < time_fit_min) ||
747  (fit_mean > time_fit_max)) {
748  B2INFO("x_bin = " << x_bin);
749  B2INFO("fit mean, default mean = " << fit_mean << ", " << default_mean);
750  B2INFO("fit mean unc, default mean unc = " << fit_mean_unc << ", " << default_mean_unc);
751  B2INFO("fit sigma, default sigma = " << fit_sigma << ", " << default_sigma);
752 
753  B2INFO("crystal fit mean - hist mean = " << meanDiff);
754  B2INFO("fit mean unc. - hist mean unc. = " << meanUncDiff);
755  B2INFO("fit sigma - hist sigma = " << sigmaDiff);
756 
757  B2INFO("fit_mean = " << fit_mean);
758  B2INFO("time_fit_min = " << time_fit_min);
759  B2INFO("time_fit_max = " << time_fit_max);
760 
761  if (fabs(meanDiff) > 10) B2INFO("fit mean diff too large");
762  if (fabs(meanUncDiff) > 10) B2INFO("fit mean unc diff too large");
763  if (fabs(sigmaDiff) > 10) B2INFO("fit mean sigma diff too large");
764  if (fit_mean_unc > 3) B2INFO("fit mean unc too large");
765  if (fit_sigma < 0.1) B2INFO("fit sigma too small");
766 
767  } else {
768  good_fit = true;
769  }
770 
771 
772  int numEntries = h_time->GetEntries();
773  /* If number of entries in histogram is greater than X then use the statistical information
774  from the data otherwise leave crystal uncalibrated. Histograms are still shown though.
775  ALSO require the that fits are good. */
776  if ((numEntries >= minNumEntries) && good_fit) {
777  clusterTime_mean = fit_mean;
778  clusterTime_mean_unc = fit_mean_unc;
779  clusterTime_sigma = fit_sigma;
780  } else {
781  clusterTime_mean = default_mean;
782  clusterTime_mean_unc = default_mean_unc;
783  clusterTime_sigma = default_sigma;
784  }
785 
786  if (numEntries < minNumEntries) B2INFO("Number of entries less than minimum");
787  if (numEntries == 0) B2INFO("Number of entries == 0");
788 
789  histfile->WriteTObject(h_time, (std::string("h_time_E_slice") + std::to_string(meanE)).c_str());
790  histfile->WriteTObject(h_timeMasked, (std::string("h_time_E_slice_masked") + std::to_string(meanE)).c_str());
791 
792  // store mean cluster time info in a separate histogram
793  clusterTimePeak_ClusterEnergy_varBin->SetBinContent(Ebin_counter, clusterTime_mean);
794  clusterTimePeak_ClusterEnergy_varBin->SetBinError(Ebin_counter, clusterTime_mean_unc);
795 
796  clusterTimePeakWidth_ClusterEnergy_varBin->SetBinContent(Ebin_counter, clusterTime_sigma);
797  clusterTimePeakWidth_ClusterEnergy_varBin->SetBinError(Ebin_counter, 0);
798 
799  Ebin_counter++;
800 
801  delete gaus;
802  }
803 
804 
805 
806  /***************************************************************************
807  For the user, print out some information about the peak cluster times.
808  It is sorted by the absolute value of the peak cluster time so that the
809  worst times are at the end.
810  ***************************************************************************/
811 
812  // Vector to store element with respective present index
813  vector< pair<double, int> > fitClusterTime__crystalIDBase0__pairs;
814 
815  // Prepare a vector of pairs containing the fitted cluster time and cell ID (base 0)
816  for (int cid = 0; cid < ECLElementNumbers::c_NCrystals; cid++) {
817  fitClusterTime__crystalIDBase0__pairs.push_back(make_pair(0.0, cid));
818  }
819 
820  // Inserting element in pair vector to keep track of crystal id.
821  for (int crys_id = cellIDLo; crys_id <= cellIDHi; crys_id++) {
822  fitClusterTime__crystalIDBase0__pairs[crys_id - 1] = make_pair(fabs(t_offsets[crys_id - 1]), crys_id - 1) ;
823  }
824 
825  // Sorting by the absolute value of the fitted cluster time for the crystal
826  sort(fitClusterTime__crystalIDBase0__pairs.begin(), fitClusterTime__crystalIDBase0__pairs.end());
827 
828 
829  // Print out the fitted peak cluster time values sorted by their absolute value
830  B2INFO("-------- List of the (fitted) peak cluster times sorted by their absolute value ----------");
831  B2INFO("------------------------------------------------------------------------------------------");
832  B2INFO("------------------------------------------------------------------------------------------");
833  B2INFO("Quoted # of clusters is before the cutting off of the distribution tails, cellID=1..ECLElementNumbers::c_NCrystals, crysID=0..8735");
834 
835  bool hasHitThresholdBadTimes = false ;
836  for (int iSortedTimes = 0; iSortedTimes < ECLElementNumbers::c_NCrystals; iSortedTimes++) {
837  int cid = fitClusterTime__crystalIDBase0__pairs[iSortedTimes].second ;
838  if (!hasHitThresholdBadTimes && fitClusterTime__crystalIDBase0__pairs[iSortedTimes].first > 2) {
839  B2INFO("======== |t_fit| > Xns threshold ======");
840  hasHitThresholdBadTimes = true;
841  }
842  //B2INFO("crystal ID = " << cid << ", peak clust t = " << t_offsets[cid] << " +- " << t_offsets_unc[cid] << ", # clusters = " << numClusterPerCrys[cid] << ", fabs(t) = " << fitClusterTime__crystalIDBase0__pairs[iSortedTimes].first );
843  B2INFO("cid = " << cid << ", peak clust t = " << t_offsets[cid] << " +- " << t_offsets_unc[cid] << " ns, # clust = " <<
844  numClusterPerCrys[cid] << ", good fit = " << crysHasGoodFit[cid] << ", good fit & stats = " << crysHasGoodFitandStats[cid]);
845  }
846 
847 
848 
849 
850 
851  // Print out just a subset that definitely don't look good even though they have good stats.
852  B2INFO("######## List of poor (fitted) peak cluster times sorted by their absolute value #########");
853  B2INFO("##########################################################################################");
854  B2INFO("##########################################################################################");
855 
856  for (int iSortedTimes = 0; iSortedTimes < ECLElementNumbers::c_NCrystals; iSortedTimes++) {
857  int cid = fitClusterTime__crystalIDBase0__pairs[iSortedTimes].second ;
858  if (fitClusterTime__crystalIDBase0__pairs[iSortedTimes].first > 2 && crysHasGoodFitandStats[cid]) {
859  B2INFO("WARNING: cid = " << cid << ", peak clust t = " << t_offsets[cid] << " +- " << t_offsets_unc[cid] << " ns, # clust = " <<
860  numClusterPerCrys[cid] << ", good fit = " << crysHasGoodFit[cid] << ", good fit & stats = " << crysHasGoodFitandStats[cid]);
861  }
862  }
863 
864 
865  B2INFO("~~~~~~~~");
866  B2INFO("Number of crystals with non-zero number of hits = " << numCrysWithNonZeroEntries);
867  B2INFO("Number of crystals with good quality fits = " << numCrysWithGoodFit);
868 
869 
870  clusterTimePeak_ClusterEnergy_varBin->ResetStats();
871  clusterTimePeakWidth_ClusterEnergy_varBin->ResetStats();
872 
873  histfile->WriteTObject(clusterTimePeak_ClusterEnergy_varBin, "clusterTimePeak_ClusterEnergy_varBin");
874  histfile->WriteTObject(clusterTimePeakWidth_ClusterEnergy_varBin, "clusterTimePeakWidth_ClusterEnergy_varBin");
875 
876 
877 
878  /* Fill histograms with the difference in the ts values from this iteration
879  and the previous values read in from the payload. */
880  B2INFO("Filling histograms for difference in crystal payload values and the pre-calibration values. These older values may be from a previous bucket or older reprocessing of the data.");
881  for (int crys_id = 1; crys_id <= ECLElementNumbers::c_NCrystals; crys_id++) {
882  double tsDiffCustomOld_ns = -999;
883  if (readPrevCrysPayload) {
884  tsDiffCustomOld_ns = (currentValuesCrys[crys_id - 1] - prevValuesCrys[crys_id - 1]) * TICKS_TO_NS;
885  B2INFO("Crystal " << crys_id << ": ts new merged - 'before 1st iter' merged = (" <<
886  currentValuesCrys[crys_id - 1] << " - " << prevValuesCrys[crys_id - 1] <<
887  ") ticks * " << TICKS_TO_NS << " ns/tick = " << tsDiffCustomOld_ns << " ns");
888 
889  }
890  tsNew_MINUS_tsCustomPrev__cid->SetBinContent(crys_id, tsDiffCustomOld_ns);
891  tsNew_MINUS_tsCustomPrev__cid->SetBinError(crys_id, 0);
892  tsNew_MINUS_tsCustomPrev__cid->ResetStats();
893 
894  tsNew_MINUS_tsCustomPrev->Fill(tsDiffCustomOld_ns);
895  tsNew_MINUS_tsCustomPrev->ResetStats();
896  }
897 
898  histfile->WriteTObject(tsNew_MINUS_tsCustomPrev__cid, "tsNew_MINUS_tsCustomPrev__cid");
899  histfile->WriteTObject(tsNew_MINUS_tsCustomPrev, "tsNew_MINUS_tsCustomPrev");
900 
901 
902  histfile->Close();
903 
904  B2INFO("Finished validations algorithm");
905  return c_OK;
906 }
void updateDBObjPtrs(const unsigned int event, const int run, const int experiment)
Updates any DBObjPtrs by calling update(event) for DBStore.
@ c_OK
Finished successfuly =0 in Python.
@ c_Failure
Failed =3 in Python.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
This class provides access to ECL channel map that is either a) Loaded from the database (see ecl/dbo...
The Class for ECL Geometry Parameters.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
void Mapping(int cid)
Mapping theta, phi Id.
int GetThetaID()
Get Theta Id.
static double getRF()
See m_rf.
int cellIDHi
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
int cellIDLo
Fit crystals with cellID0 in the inclusive range [cellIDLo,cellIDHi].
double meanCleanRebinFactor
Rebinning factor for mean calculation.
double clusterTimesFractionWindow_maxtime
Maximum time for window to calculate cluster time fraction, in ns.
double meanCleanCutMinFactor
After rebinning, create a mask for bins that have values less than meanCleanCutMinFactor times the ma...
bool readPrevCrysPayload
Read the previous crystal payload values for comparison.
std::string debugFilenameBase
Name of file with debug output, eclTValidationAlgorithm.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.

◆ checkPyExpRun()

bool checkPyExpRun ( PyObject *  pyObj)
inherited

Checks that a PyObject can be successfully converted to an ExpRun type.

Checks if the PyObject can be converted to ExpRun.

Definition at line 28 of file CalibrationAlgorithm.cc.

◆ convertPyExpRun()

ExpRun convertPyExpRun ( PyObject *  pyObj)
inherited

Performs the conversion of PyObject to ExpRun.

Converts the PyObject to an ExpRun. We've preoviously checked the object so this assumes a lot about the PyObject.

Definition at line 70 of file CalibrationAlgorithm.cc.

◆ execute()

CalibrationAlgorithm::EResult execute ( std::vector< Calibration::ExpRun >  runs = {},
int  iteration = 0,
IntervalOfValidity  iov = IntervalOfValidity() 
)
inherited

Runs calibration over vector of runs for a given iteration.

You can also specify the IoV to save the database payload as. By default the Algorithm will create an IoV from your requested ExpRuns, or from the overall ExpRuns of the input data if you haven't specified ExpRuns in this function.

No checks are performed to make sure that a IoV you specify matches the data you ran over, it simply labels the IoV to commit to the database later.

Definition at line 114 of file CalibrationAlgorithm.cc.

◆ getCollectorName()

std::string getCollectorName ( ) const
inlineinherited

Alias for prefix.

For convenience and less writing, we say developers to set this to default collector module name in constructor of base class. One can however use the dublets of collector+algorithm multiple times with different settings. To bind these together correctly, the prefix has to be set the same for algo and collector. So we call the setter setPrefix rather than setModuleName or whatever. This getter will work out of the box for default cases -> return the name of module you have to add to your path to collect data for this algorihtm.

Definition at line 164 of file CalibrationAlgorithm.h.

◆ setInputFileNames() [1/2]

void setInputFileNames ( PyObject *  inputFileNames)
inherited

Set the input file names used for this algorithm from a Python list.

Set the input file names used for this algorithm and resolve the wildcards.

Definition at line 166 of file CalibrationAlgorithm.cc.

◆ setInputFileNames() [2/2]

void setInputFileNames ( std::vector< std::string >  inputFileNames)
protectedinherited

Set the input file names used for this algorithm.

Set the input file names used for this algorithm and resolve the wildcards.

Definition at line 194 of file CalibrationAlgorithm.cc.

Member Data Documentation

◆ meanCleanCutMinFactor

double meanCleanCutMinFactor

After rebinning, create a mask for bins that have values less than meanCleanCutMinFactor times the maximum bin value.

Expand mask and apply to non-rebinned histogram.

Definition at line 49 of file eclTValidationAlgorithm.h.


The documentation for this class was generated from the following files: