Belle II Software  release-06-02-00
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
 

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)
 

Detailed Description

Validate the ecl timing calibrations using a hadronic event selection.

Definition at line 27 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 82 of file eclTValidationAlgorithm.cc.

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

◆ 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()

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.

◆ 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 45 of file eclTValidationAlgorithm.h.


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