Belle II Software development
eclNOptimalAlgorithm Class Reference

Algorithm that works with eclNOptimalCollector to find the number of crystals to be summed to get the best energy resolution for each test energy and for each group of crystals (8 groups per thetaID in the barrel). More...

#include <eclNOptimalAlgorithm.h>

Inheritance diagram for eclNOptimalAlgorithm:
CalibrationAlgorithm

Public Types

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

Public Member Functions

 eclNOptimalAlgorithm ()
 ..Constructor
 
virtual ~eclNOptimalAlgorithm ()
 ..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.
 
Calibration::ExpRun convertPyExpRun (PyObject *pyObj)
 Performs the conversion of PyObject to ExpRun.
 
std::string getCollectorName () const
 Alias for prefix.
 
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.
 
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.
 
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 algorithm (set by developers in constructor)
 
bool loadInputJson (const std::string &jsonString)
 Load the m_inputJson variable from a string (useful from Python interface). The return 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>.
 

Protected Member Functions

virtual EResult calibrate () override
 ..Run algorithm on events
 
void setInputFileNames (std::vector< std::string > inputFileNames)
 Set the input file names used for this algorithm.
 
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 internally 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

Algorithm that works with eclNOptimalCollector to find the number of crystals to be summed to get the best energy resolution for each test energy and for each group of crystals (8 groups per thetaID in the barrel).

Also finds the corresponding energy bias from beam backgrounds, and the peak fraction of energy contained in the crystals.

Contact Chris Hearty heart.nosp@m.y@ph.nosp@m.ysics.nosp@m..ubc.nosp@m..ca for questions or concerns Find optimal number of crystals to sum for cluster energy

Definition at line 30 of file eclNOptimalAlgorithm.h.

Member Enumeration Documentation

◆ EResult

enum EResult
inherited

The result of calibration.

Enumerator
c_OK 

Finished successfully =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.

40 {
41 c_OK,
42 c_Iterate,
44 c_Failure,
46 };
@ c_OK
Finished successfully =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.

Constructor & Destructor Documentation

◆ eclNOptimalAlgorithm()

..Constructor

Definition at line 210 of file eclNOptimalAlgorithm.cc.

210 : CalibrationAlgorithm("eclNOptimalCollector")
211{
213 "Use single gamma MC to find the optimal number of crystals to sum for cluster energy"
214 );
215}
Base class for calibration algorithms.
void setDescription(const std::string &description)
Set algorithm description (in constructor)

◆ ~eclNOptimalAlgorithm()

virtual ~eclNOptimalAlgorithm ( )
inlinevirtual

..Destructor

Definition at line 37 of file eclNOptimalAlgorithm.h.

37{}

Member Function Documentation

◆ boundaryFindingSetup()

virtual void boundaryFindingSetup ( std::vector< Calibration::ExpRun >  ,
int   
)
inlineprotectedvirtualinherited

If you need to make some changes to your algorithm class before 'findPayloadBoundaries' is run, make them in this function.

Reimplemented in TestBoundarySettingAlgorithm, TestCalibrationAlgorithm, PXDAnalyticGainCalibrationAlgorithm, PXDValidationAlgorithm, SVD3SampleCoGTimeCalibrationAlgorithm, SVD3SampleELSTimeCalibrationAlgorithm, and SVDCoGTimeCalibrationAlgorithm.

Definition at line 252 of file CalibrationAlgorithm.h.

252{};

◆ boundaryFindingTearDown()

virtual void boundaryFindingTearDown ( )
inlineprotectedvirtualinherited

Put your algorithm back into a state ready for normal execution if you need to.

Definition at line 257 of file CalibrationAlgorithm.h.

257{};

◆ calibrate()

CalibrationAlgorithm::EResult calibrate ( )
overrideprotectedvirtual

..Run algorithm on events

Implements CalibrationAlgorithm.

Definition at line 219 of file eclNOptimalAlgorithm.cc.

220{
221
222 //-----------------------------------------------------------------------------------
223 //-----------------------------------------------------------------------------------
224 //..Algorithm set up
225
226 //..Get started
227 B2INFO("Starting eclNOptimalAlgorithm");
228 gROOT->SetBatch();
229
230 //-----------------------------------------------------------------------------------
231 //..Read in parameter histograms created by the collector and fix normalization
232 auto inputParameters = getObjectPtr<TH1F>("inputParameters");
233 const int lastBin = inputParameters->GetNbinsX();
234 const double scale = inputParameters->GetBinContent(lastBin); // number of times inputParameters was filled
235 for (int ib = 1; ib < lastBin; ib++) {
236 double param = inputParameters->GetBinContent(ib);
237 inputParameters->SetBinContent(ib, param / scale);
238 inputParameters->SetBinError(ib, 0.);
239 }
240
241 //..Also read in and normalize group number for each cellID
242 auto groupNumberOfEachCellID = getObjectPtr<TH1F>("groupNumberOfEachCellID");
243 for (int ic = 1; ic <= ECLElementNumbers::c_NCrystals; ic++) {
244 const double groupNum = groupNumberOfEachCellID->GetBinContent(ic);
245 groupNumberOfEachCellID->SetBinContent(ic, groupNum / scale);
246 groupNumberOfEachCellID->SetBinError(ic, 0.);
247 }
248
249 //..Couple of diagnostic histograms
250 auto entriesPerThetaIdEnergy = getObjectPtr<TH2F>("entriesPerThetaIdEnergy");
251 auto mcEnergyDiff = getObjectPtr<TH2F>("mcEnergyDiff");
252
253 //..Write these to disk.
254 TFile* histFile = new TFile("eclNOptimalAlgorithm.root", "recreate");
255 inputParameters->Write();
256 groupNumberOfEachCellID->Write();
257 entriesPerThetaIdEnergy->Write();
258 mcEnergyDiff->Write();
259
260 //-----------------------------------------------------------------------------------
261 //..Parameters from the inputParameters histogram
262 const int nCrystalGroups = (int)(inputParameters->GetBinContent(1) + 0.5);
263 const int nEnergies = (int)(inputParameters->GetBinContent(2) + 0.5);
264 const int nLeakReg = 3; // 3 regions forward barrel backward
265 auto generatedE = create2Dvector<float>(nLeakReg, nEnergies);
266 int bin = 2; // bin 1 = nCrystalGroups, bin 2 = nEnergies, bin 3 = first energy
267 for (int ireg = 0; ireg < nLeakReg; ireg++) {
268 for (int ie = 0; ie < nEnergies; ie++) {
269 bin++;
270 generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
271 }
272 }
273 B2INFO("nCrystalGroups = " << nCrystalGroups);
274
275 //-----------------------------------------------------------------------------------
276 //..Experiment and run number, for reading existing payload from database
277 const auto exprun = getRunList();
278 const int iExp = exprun[0].first;
279 const int iRun = exprun[0].second;
280 B2INFO("Experiment, run = " << iExp << ", " << iRun);
281
283 // simulate the initialize() phase where we can register objects in the DataStore
285 evtPtr.registerInDataStore();
287 // now construct the event metadata
288 evtPtr.construct(1, iRun, iExp);
289 // and update the database contents
290 DBStore& dbstore = DBStore::Instance();
291 dbstore.update();
292 // this is only needed it the payload might be intra-run dependent,
293 // that is if it might change during one run as well
294 dbstore.updateEvent();
295
296 //-----------------------------------------------------------------------------------
297 //..Get existing payload
298 B2INFO("eclNOptimalAlgorithm: reading previous payload");
299 DBObjPtr<ECLnOptimal> existingECLnOptimal;
300 std::vector<float> eUpperBoundariesFwdPrev = existingECLnOptimal->getUpperBoundariesFwd();
301 std::vector<float> eUpperBoundariesBrlPrev = existingECLnOptimal->getUpperBoundariesBrl();
302 std::vector<float> eUpperBoundariesBwdPrev = existingECLnOptimal->getUpperBoundariesBwd();
303 TH2F nOptimal2DPrev = existingECLnOptimal->getNOptimal();
304 nOptimal2DPrev.SetName("nOptimal2DPrev");
305
306 const int nPrevE = eUpperBoundariesFwdPrev.size();
307 B2INFO("Energy upper boundaries from previous payload for each region");
308 for (int ie = 0; ie < nPrevE; ie++) {
309 B2INFO(" " << ie << " " << eUpperBoundariesFwdPrev[ie] << " " << eUpperBoundariesBrlPrev[ie] << " " << eUpperBoundariesBwdPrev[ie]);
310 }
311
312 //..Write out
313 histFile->cd();
314 nOptimal2DPrev.Write();
315
316 //-----------------------------------------------------------------------------------
317 //..Use existing payload to get the initial nOptimal value for each group
318 std::vector< std::vector<int> > initialNOptimal;
319 initialNOptimal.resize(nCrystalGroups, std::vector<int>(nEnergies, 0));
320
321 //..Couple of histograms of relevant quantities
322 TH2F* nInitialPerGroup = new TH2F("nInitialPerGroup", "initial nCrys, energy point vs group number;group;energy point",
323 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
324 TH2F* generatedEPerGroup = new TH2F("generatedEPerGroup", "Generated energy, energy point vs group number;group;energy point",
325 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
326 TH1F* firstCellIDPerGroup = new TH1F("firstCellIDPerGroup", "First cellID of group;group", nCrystalGroups, 0., nCrystalGroups);
327 TH1F* nCrystalsPerGroup = new TH1F("nCrystalsPerGroup", "Number of crystals per group;group", nCrystalGroups, 0., nCrystalGroups);
328 std::vector<int> crystalsPerGroup;
329 crystalsPerGroup.resize(nCrystalGroups, 0);
330
331 const int nCrystalsTotal = ECLElementNumbers::c_NCrystals;
332 for (int cellID = 1; cellID <= nCrystalsTotal; cellID++) {
333
334 //..Group number of this cellID
335 int iGroup = (int)(0.5 + groupNumberOfEachCellID->GetBinContent(cellID));
336
337 //..Energy boundaries of previous payload, which depend on the ECL region
338 std::vector<float> eUpperBoundariesPrev = eUpperBoundariesBrlPrev;
339 int iRegion = 1; // barrel
340 if (ECLElementNumbers::isForward(cellID)) {
341 eUpperBoundariesPrev = eUpperBoundariesFwdPrev;
342 iRegion = 0;
343 } else if (ECLElementNumbers::isBackward(cellID)) {
344 eUpperBoundariesPrev = eUpperBoundariesBwdPrev;
345 iRegion = 2;
346 }
347
348 //..For each test energy, get the nOptimal for this cellID from previous payload.
349 // iEnergy is the energy bin of the previous payload, which is equal to ie if
350 // the test energies have not changed.
351 for (int ie = 0; ie < nEnergies; ie++) {
352 float energy = generatedE[iRegion][ie];
353 int iEnergy = 0;
354 while (energy > eUpperBoundariesPrev[iEnergy]) {iEnergy++;}
355 initialNOptimal[iGroup][ie] = (int)(0.5 + nOptimal2DPrev.GetBinContent(iGroup + 1, iEnergy + 1));
356 if (ie != iEnergy) {
357 B2INFO("ie iEnergy mismatch: cellID " << cellID << " ie " << ie << " iEnergy " << iEnergy);
358 }
359
360 //..Store in diagnostic histograms
361 generatedEPerGroup->SetBinContent(iGroup + 1, ie + 1, energy);
362 nInitialPerGroup->SetBinContent(iGroup + 1, ie + 1, initialNOptimal[iGroup][ie]);
363 }
364
365 //..Count crystals in this group
366 if (crystalsPerGroup[iGroup] == 0) {
367 firstCellIDPerGroup->SetBinContent(iGroup + 1, cellID);
368 firstCellIDPerGroup->SetBinError(iGroup + 1, 0.);
369 }
370 crystalsPerGroup[iGroup]++;
371 }
372
373 //..Write out these histograms
374 for (int ig = 0; ig < nCrystalGroups; ig++) {
375 nCrystalsPerGroup->SetBinContent(ig + 1, crystalsPerGroup[ig]);
376 nCrystalsPerGroup->SetBinError(ig + 1, 0.);
377 }
378 histFile->cd();
379 nInitialPerGroup->Write();
380 generatedEPerGroup->Write();
381 firstCellIDPerGroup->Write();
382 nCrystalsPerGroup->Write();
383 B2INFO("Wrote initial diagnostic histograms");
384
385 //-----------------------------------------------------------------------------------
386 //..Write out all the 2D histograms
387 for (int ig = 0; ig < nCrystalGroups; ig++) {
388 for (int ie = 0; ie < nEnergies; ie++) {
389 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
390 auto eSum = getObjectPtr<TH2F>(name);
391 name = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
392 auto biasSum = getObjectPtr<TH2F>(name);
393 histFile->cd();
394 eSum->Write();
395 biasSum->Write();
396 }
397 }
398 B2INFO("Wrote eSum and biasSum histograms");
399
400 //-----------------------------------------------------------------------------------
401 //-----------------------------------------------------------------------------------
402 //..Find nOptimal and corresponding bias for each group/energy point
403
404 //..Histograms to store fit results
405 TH2F* nOptimalPerGroup = new TH2F("nOptimalPerGroup", "nOptimal;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies,
406 0., nEnergies);
407 TH2F* changeInNOptimal = new TH2F("changeInNOptimal", "nOptimal minus nInitial;group;energy point", nCrystalGroups, 0.,
408 nCrystalGroups, nEnergies, 0., nEnergies);
409 TH2F* binsInFit = new TH2F("binsInFit", "Number of bins used in Novo fit, energy point vs group number;group;energy point",
410 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
411 TH2F* maxEntriesPerHist = new TH2F("maxEntriesPerHist",
412 "Max entries in histogram bin, energy point vs group number;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
413 nEnergies);
414 TH2F* peakPerGroup = new TH2F("peakPerGroup", "peak energy, energy point vs group number;group;energy point", nCrystalGroups, 0.,
415 nCrystalGroups, nEnergies, 0., nEnergies);
416 TH2F* effSigmaPerGroup = new TH2F("effSigmaPerGroup", "fit effective sigma, energy point vs group number;group;energy point",
417 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
418 TH2F* etaPerGroup = new TH2F("etaPerGroup", "fit eta, energy point vs group number;group;energy point", nCrystalGroups, 0.,
419 nCrystalGroups, nEnergies, 0., nEnergies);
420 TH2F* fitProbPerGroup = new TH2F("fitProbPerGroup", "fit probability, energy point vs group number;group;energy point",
421 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
422 TH2F* resolutionPerGroup = new TH2F("resolutionPerGroup", "resolution/(peak-bias), energy point vs group number;group;energy point",
423 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
424 TH2F* fracContainedEPerGroup = new TH2F("fracContainedEPerGroup",
425 "peak fraction of energy contained in nOpt crystals;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
426 nEnergies);
427
428 TH2F* biasPerGroup = new TH2F("biasPerGroup", "bias (GeV), energy point vs group number;group;energy point", nCrystalGroups, 0.,
429 nCrystalGroups, nEnergies, 0., nEnergies);
430 TH2F* sigmaBiasPerGroup = new TH2F("sigmaBiasPerGroup", "sigma bias (GeV), energy point vs group number;group;energy point",
431 nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
432
433 TH2F* existingResolutionPerGroup = new TH2F("existingResolutionPerGroup",
434 "existing resolution/peak, energy point vs group number;group;energy point", nCrystalGroups, 0., nCrystalGroups, nEnergies, 0.,
435 nEnergies);
436
437 //-----------------------------------------------------------------------------------
438 //..Loop over all groups and energy points
439
440 //..Fit function Novosibirsk (xmin, xmax, nparameters)
441 TF1* func = new TF1("eclNOptimalNovo", eclNOptimalNovo, 0.4, 1.4, 4);
442 func->SetLineColor(kRed);
443
444 for (int ig = 0; ig < nCrystalGroups; ig++) {
445 for (int ie = 0; ie < nEnergies; ie++) {
446
447 //..Read in the energy sum vs nCrys for the group and energy point
448 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ie);
449 auto eSum = getObjectPtr<TH2F>(name);
450
451 //..Also the corresponding bias histogram
452 std::string biasName = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ie);
453 auto biasSum = getObjectPtr<TH2F>(biasName);
454
455 //..And the generated energy
456 const double genEnergy = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
457
458 //-----------------------------------------------------------------------------------
459 //..Extract the eSum distributions for various nCrys (=nCrysSumToFit) and fit them
460 // to find the best resolution.
461 const int nCrysBins = eSum->GetNbinsX();
462 const int nCrysMax = nCrysBins - 1; // the last bin is eSum for existing reco
463
464 //..Items to be found and stored in the following while loop
465 TH1D* eSumOpt{nullptr}; // 1D projection with the best resolution
466 int nFitBins = 0; // keep track of how many bins are included in the fit
467 double bestFractionalResolution = 999.;
468 double existingFractionalResolution = 999.;
469 double biasForNOpt = 0.; // bias from backgrounds for optimal nCrys
470 double biasSigmaForNOpt = 0.; // sigma bias
471 int nOpt = 0;
472
473 //..Start with the previous optimal value, then try fewer, if possible, then
474 // try more, if possible. End by finding the resolution for the
475 // reconstruction used to produce the single photon MC samples.
476 const int initialnCrysSumToFit = initialNOptimal[ig][ie];
477 int nCrysSumToFit = initialnCrysSumToFit;
478 while (nCrysSumToFit > 0) {
479
480 TH1D* hEnergy = (TH1D*)eSum->ProjectionY("hEnergy", nCrysSumToFit, nCrysSumToFit);
481 TString newName = name + "_" + std::to_string(nCrysSumToFit);
482 hEnergy->SetName(newName);
483
484 //------------------------------------------------------------------------------
485 //..Check stats, and rebin as required
486 const double minESumEntries = 800.;
487 if (hEnergy->GetEntries() < minESumEntries) {
488 B2INFO("Insufficient entries in eSum: " << hEnergy->GetEntries() << " for group " << ig << " energy point " << ie);
489 histFile->Close();
490 B2INFO("closed histFile; quitting");
491 return c_NotEnoughData;
492 }
493
494 //..Rebin if the maximum bin has too few entries
495 const double minPeakValue = 300.;
496 while (hEnergy->GetMaximum() < minPeakValue) {hEnergy->Rebin(2);}
497
498 //------------------------------------------------------------------------------
499 //..Find 50% range for Novo fit to find peak value. Rebin if there are too many bins
500 double fitFraction = 0.5;
501 int iLo = 1;
502 int iHi = hEnergy->GetNbinsX();
503 const int maxBinsToFit = 59;
504 while ((iHi - iLo) > maxBinsToFit) {
505 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
506 iLo = iBins[0];
507 iHi = iBins[1];
508 if ((iHi - iLo) > maxBinsToFit) {hEnergy->Rebin(2);}
509 }
510
511 //..Fit range
512 const int nBinsX = hEnergy->GetNbinsX();
513 const double fitlow = hEnergy->GetBinLowEdge(iLo);
514 const double fithigh = hEnergy->GetBinLowEdge(iHi + 1);
515 double sum = hEnergy->Integral(iLo, iHi);
516 B2INFO("group " << ig << " ie " << ie << " name " << name << " newName " << newName << " nBinsX " << nBinsX << " 50% target: " <<
517 fitFraction * hEnergy->GetEntries() << " iLo: " << iLo << " iHi: " << iHi << " sum: " << sum);
518
519 //------------------------------------------------------------------------------
520 //..Set up the fit. Initial values
521 double normalization = hEnergy->GetMaximum();
522 const int iMax = hEnergy->GetMaximumBin();
523 double peak = hEnergy->GetBinLowEdge(iMax);
524 double effSigma = 0.5 * (fithigh - fitlow);
525 double eta = 0.4; // typical value for good fits
526
527 func->SetParameters(normalization, peak, effSigma, eta);
528 func->SetParNames("normalization", "peak", "effSigma", "eta");
529
530 //..Parameter ranges
531 func->SetParLimits(1, fitlow, fithigh);
532 func->SetParLimits(2, 0., 2.*effSigma);
533 func->SetParLimits(3, -1, 3);
534
535 //..If there is only one crystal, fix eta to a special value
536 if (nCrysSumToFit == 1) {
537 const double etaForOneCrystal = 0.95;
538 func->SetParameter(3, etaForOneCrystal);
539 func->SetParLimits(3, etaForOneCrystal, etaForOneCrystal);
540 }
541
542 //..Fit
543 hEnergy->Fit(func, "LIEQ", "", fitlow, fithigh);
544
545 //------------------------------------------------------------------------------
546 //..Find the bias = sum of CalDigits minus sum of MC truth for this nCrys
547 // Bias is the mid point of the minimum range that contains 68.3% of events
548 double bias = 0.;
549 double biasSigma = 0.;
550 if (nCrysSumToFit < nCrysBins) {
551 TH1D* hBias = (TH1D*)biasSum->ProjectionY("hBias", nCrysSumToFit, nCrysSumToFit);
552 fitFraction = 0.683;
553 std::vector<int> jBins = eclNOptimalFitRange(hBias, fitFraction);
554 const double lowEdge = hBias->GetBinLowEdge(jBins[0]);
555 const double highEdge = hBias->GetBinLowEdge(jBins[1] + 1);
556 bias = 0.5 * (lowEdge + highEdge);
557 biasSigma = 0.5 * (highEdge - lowEdge);
558 }
559
560 //------------------------------------------------------------------------------
561 //..Resolution is half the smallest range that contains 68.3% of events.
562 // Speed up the process by only looking among the range that contains 75%.
563 fitFraction = 0.75;
564 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
565 const double resolution = eclNOptimalResolution(hEnergy, iBins[0], iBins[1]);
566
567 //..Correct for the bias from background when calculating fractional resolution.
568 // Note that variable resolution = (E_hi - E_lo)/E_gen is dimensionless,
569 // variable peak = E_peak/E_gen is dimensionless, but variable bias has units GeV
570 peak = func->GetParameter(1);
571 const double fractionalResolution = resolution / (peak - bias / genEnergy);
572
573 //..Record information if this is the best resolution (and is not the existing
574 // resolution case).
575 if (fractionalResolution < bestFractionalResolution and nCrysSumToFit != nCrysBins) {
576 bestFractionalResolution = fractionalResolution;
577 nOpt = nCrysSumToFit;
578 eSumOpt = hEnergy;
579 nFitBins = 1 + iHi - iLo;
580 biasForNOpt = bias;
581 biasSigmaForNOpt = biasSigma;
582 } else if (nCrysSumToFit == nCrysBins) {
583 existingFractionalResolution = fractionalResolution;
584 }
585
586 //------------------------------------------------------------------------------
587 //..Logic to decide what to do next
588
589 //..keep checking different N until resolution is this much worse than best value
590 const double resTolerance = 1.05;
591 if (nCrysSumToFit == initialnCrysSumToFit) {
592
593 //..After testing the previous nOptimal, try one fewer if possible
594 if (nCrysSumToFit > 1) {
595 nCrysSumToFit--;
596 } else {
597 nCrysSumToFit++;
598 }
599 } else if (nCrysSumToFit < initialnCrysSumToFit) {
600
601 //..Trying fewer crystals. If this is the best so far, try one
602 // fewer, if possible. Otherwise, try more crystals.
603 if (fractionalResolution < resTolerance* bestFractionalResolution and nCrysSumToFit > 1) {
604 nCrysSumToFit--;
605 } else if (initialnCrysSumToFit != nCrysMax) {
606 nCrysSumToFit = initialnCrysSumToFit + 1;
607 } else {
608 nCrysSumToFit = nCrysBins;
609 }
610 } else if (nCrysSumToFit < nCrysBins) { // nCrysSumToFit is always > initialnCrysSumToFit
611
612 //..Trying more crystals. If this is the best, try one more, if
613 // possible. Otherwise, do the current reconstruction case (nCrysSumToFit = nCrysBins)
614 if (fractionalResolution < resTolerance * bestFractionalResolution and nCrysSumToFit < nCrysMax) {
615 nCrysSumToFit++;
616 } else {
617 nCrysSumToFit = nCrysBins;
618 }
619 } else if (nCrysSumToFit == nCrysBins) {
620
621 //..This is the current resolution case, so we are done
622 nCrysSumToFit = 0;
623 }
624
625 } // while(nCrysSumToFit > 0)
626 B2INFO(" ig: " << ig << " ie: " << ie << " initial nOpt: " << initialnCrysSumToFit << " final nOpt: " << nOpt);
627
628 //-----------------------------------------------------------------------------------
629 //..Store everything in diagnostic histograms
630
631 //..Extract the function from the nOptimal histogram
632 TF1* funcOpt = (TF1*)eSumOpt->GetFunction("eclNOptimalNovo");
633
634 nOptimalPerGroup->SetBinContent(ig + 1, ie + 1, nOpt);
635
636 changeInNOptimal->SetBinContent(ig + 1, ie + 1, nOpt - initialnCrysSumToFit);
637 changeInNOptimal->SetBinError(ig + 1, ie + 1, 0.);
638
639 biasPerGroup->SetBinContent(ig + 1, ie + 1, biasForNOpt);
640 biasPerGroup->SetBinError(ig + 1, ie + 1, 0.);
641
642 sigmaBiasPerGroup->SetBinContent(ig + 1, ie + 1, biasSigmaForNOpt);
643 sigmaBiasPerGroup->SetBinError(ig + 1, ie + 1, 0.);
644
645 binsInFit->SetBinContent(ig + 1, ie + 1, nFitBins);
646 binsInFit->SetBinError(ig + 1, ie + 1, 0.);
647
648 maxEntriesPerHist->SetBinContent(ig + 1, ie + 1, eSumOpt->GetMaximum());
649 maxEntriesPerHist->SetBinError(ig + 1, ie + 1, 0.);
650
651 const double peakOpt = funcOpt->GetParameter(1);
652 const double peakOptUnc = funcOpt->GetParError(1);
653 peakPerGroup->SetBinContent(ig + 1, ie + 1, peakOpt);
654 peakPerGroup->SetBinError(ig + 1, ie + 1, peakOptUnc);
655
656 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
657 const double peakCor = peakOpt - biasForNOpt / genE;
658 fracContainedEPerGroup->SetBinContent(ig + 1, ie + 1, peakCor);
659 fracContainedEPerGroup->SetBinError(ig + 1, ie + 1, peakOptUnc);
660
661 effSigmaPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetParameter(2));
662 effSigmaPerGroup->SetBinError(ig + 1, ie + 1, funcOpt->GetParError(2));
663
664 etaPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetParameter(3));
665 etaPerGroup->SetBinError(ig + 1, ie + 1, funcOpt->GetParError(3));
666
667 fitProbPerGroup->SetBinContent(ig + 1, ie + 1, funcOpt->GetProb());
668 fitProbPerGroup->SetBinError(ig + 1, ie + 1, 0.);
669
670 resolutionPerGroup->SetBinContent(ig + 1, ie + 1, bestFractionalResolution);
671 resolutionPerGroup->SetBinError(ig + 1, ie + 1, 0.);
672
673 existingResolutionPerGroup->SetBinContent(ig + 1, ie + 1, existingFractionalResolution);
674 existingResolutionPerGroup->SetBinError(ig + 1, ie + 1, 0.);
675
676 //..Write out to disk
677 histFile->cd();
678 eSumOpt->Write();
679
680 }
681 }
682
683 //-----------------------------------------------------------------------------------
684 //..Write out fit summary histograms
685 histFile->cd();
686 changeInNOptimal->Write();
687 biasPerGroup->Write();
688 sigmaBiasPerGroup->Write();
689 binsInFit->Write();
690 maxEntriesPerHist->Write();
691 peakPerGroup->Write();
692 effSigmaPerGroup->Write();
693 etaPerGroup->Write();
694 fitProbPerGroup->Write();
695 resolutionPerGroup->Write();
696 existingResolutionPerGroup->Write();
697 fracContainedEPerGroup->Write();
698
699 B2INFO("Wrote fit summary histograms");
700
701 //-----------------------------------------------------------------------------------
702 //-----------------------------------------------------------------------------------
703 //..Find bias and fraction contained energy in nOpt crystals for the energy bins
704 // adjacent to each of the group / energy points.
705
706 //..Histograms to store results
707 TH2F* fracContainedAdjacent[2];
708 TH2F* biasAdjacent[2];
709 const TString updown[2] = {"lower", "higher"};
710 for (int ia = 0; ia < 2; ia++) {
711 TString name = "fracContainedAdjacent_" + std::to_string(ia);
712 TString title = "peak fraction of energy contained in nOpt crystals, " + updown[ia] + " E;group;energy point";
713 fracContainedAdjacent[ia] = new TH2F(name, title, nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
714
715 name = "biasAdjacent_" + std::to_string(ia);
716 title = "bias (GeV) in nOpt crystals, " + updown[ia] + " E;group;energy point";
717 biasAdjacent[ia] = new TH2F(name, title, nCrystalGroups, 0., nCrystalGroups, nEnergies, 0., nEnergies);
718 }
719
720 //-----------------------------------------------------------------------------------
721 //..Loop over all groups and energy points
722 for (int ig = 0; ig < nCrystalGroups; ig++) {
723 for (int ie = 0; ie < nEnergies; ie++) {
724
725 //..nOpt for this point
726 const int nOpt = nOptimalPerGroup->GetBinContent(ig + 1, ie + 1);
727
728 //-----------------------------------------------------------------------------------
729 //..Two adjacent energy points
730 for (int ia = 0; ia < 2; ia++) {
731 const int ieAdj = ie + 2 * ia - 1; // ie +/- 1
732 double eFracPeakAdj = fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1);
733 double biasAdj = biasPerGroup->GetBinContent(ig + 1, ie + 1);
734
735 //..Need this to be a valid energy point to do anything more
736 if (ieAdj >= 0 and ieAdj < nEnergies) {
737
738 //----------------------------------------------------------------------------
739 //..Find the bias
740 std::string biasName = "biasSum_" + std::to_string(ig) + "_" + std::to_string(ieAdj);
741 auto biasSum = getObjectPtr<TH2F>(biasName);
742 TH1D* hBias = (TH1D*)biasSum->ProjectionY("hBias", nOpt, nOpt);
743
744 //..Bias is the mid-point of the range containing 68.3% of events
745 double fitFraction = 0.683;
746 std::vector<int> jBins = eclNOptimalFitRange(hBias, fitFraction);
747 const double lowEdge = hBias->GetBinLowEdge(jBins[0]);
748 const double highEdge = hBias->GetBinLowEdge(jBins[1] + 1);
749 biasAdj = 0.5 * (lowEdge + highEdge);
750
751 //----------------------------------------------------------------------------
752 //..Get the eSum distribution to be fit
753 std::string name = "eSum_" + std::to_string(ig) + "_" + std::to_string(ieAdj);
754 auto eSum = getObjectPtr<TH2F>(name);
755 TH1D* hEnergy = (TH1D*)eSum->ProjectionY("hEnergy", nOpt, nOpt);
756 TString newName = name + "_" + std::to_string(nOpt);
757 hEnergy->SetName(newName);
758
759 B2INFO("fit adj ig = " << ig << " ie = " << ie << " ia = " << ia << " " << newName << " " << " entries = " << hEnergy->GetEntries()
760 << " GetMaxiumum = " << hEnergy->GetMaximum());
761
762 //..Rebin if the maximum bin has too few entries
763 const double minPeakValue = 300.;
764 while (hEnergy->GetMaximum() < minPeakValue) {hEnergy->Rebin(2);}
765
766 //..Find 50% range for Novo fit. Rebin if there are too many bins
767 fitFraction = 0.5;
768 int iLo = 1;
769 int iHi = hEnergy->GetNbinsX();
770 const int maxBinsToFit = 59;
771 while ((iHi - iLo) > maxBinsToFit) {
772 std::vector<int> iBins = eclNOptimalFitRange(hEnergy, fitFraction);
773 iLo = iBins[0];
774 iHi = iBins[1];
775 if ((iHi - iLo) > maxBinsToFit) {hEnergy->Rebin(2);}
776 }
777
778 //..Fit parameters
779 const double fitlow = hEnergy->GetBinLowEdge(iLo);
780 const double fithigh = hEnergy->GetBinLowEdge(iHi + 1);
781 double normalization = hEnergy->GetMaximum();
782 const int iMax = hEnergy->GetMaximumBin();
783 double peak = hEnergy->GetBinLowEdge(iMax);
784 double effSigma = 0.5 * (fithigh - fitlow);
785 double eta = 0.4; // typical value for good fits
786
787 func->SetParameters(normalization, peak, effSigma, eta);
788 func->SetParNames("normalization", "peak", "effSigma", "eta");
789
790 //..Parameter ranges
791 func->SetParLimits(1, fitlow, fithigh);
792 func->SetParLimits(2, 0., 2.*effSigma);
793 func->SetParLimits(3, -1, 3);
794
795 //..If there is only one crystal, fix eta to a special value
796 if (nOpt == 1) {
797 const double etaForOneCrystal = 0.95;
798 func->SetParameter(3, etaForOneCrystal);
799 func->SetParLimits(3, etaForOneCrystal, etaForOneCrystal);
800 }
801
802 //..Fit
803 hEnergy->Fit(func, "LIEQ", "", fitlow, fithigh);
804
805 //..Fraction of contained energy = peak of fit corrected for bias
806 const double peakAdj = func->GetParameter(1);
807 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ieAdj + 1);
808 eFracPeakAdj = peakAdj - biasAdj / genE;
809 }
810
811 //..Store
812 fracContainedAdjacent[ia]->SetBinContent(ig + 1, ie + 1, eFracPeakAdj);
813 biasAdjacent[ia]->SetBinContent(ig + 1, ie + 1, biasAdj);
814 } // loop over adjacent energy points
815
816 B2INFO(" ig " << ig << " ie " << ie << " eFrac (3): "
817 << fracContainedAdjacent[0]->GetBinContent(ig + 1, ie + 1) << " "
818 << fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1) << " "
819 << fracContainedAdjacent[1]->GetBinContent(ig + 1, ie + 1)
820 << " biases (3): "
821 << biasAdjacent[0]->GetBinContent(ig + 1, ie + 1) << " "
822 << biasPerGroup->GetBinContent(ig + 1, ie + 1) << " "
823 << biasAdjacent[1]->GetBinContent(ig + 1, ie + 1)
824 );
825 } // loop over ie
826 } // loop over ig
827
828 //-----------------------------------------------------------------------------------
829 //..Write out histograms with results for adjacent energy points
830 histFile->cd();
831 fracContainedAdjacent[0]->Write();
832 fracContainedAdjacent[1]->Write();
833 biasAdjacent[0]->Write();
834 biasAdjacent[1]->Write();
835
836 B2INFO("Wrote adjacent energy point summary histograms");
837
838
839 //-----------------------------------------------------------------------------------
840 //-----------------------------------------------------------------------------------
841 //..Prepare the payload contents
842 //..Upper energy boundaries are the mid-point of the log energies
843 auto boundaryE = create2Dvector<float>(nLeakReg, nEnergies);
844 for (int ireg = 0; ireg < nLeakReg; ireg++) {
845 B2INFO("Generated energies and boundaries for region = " << ireg);
846 for (int ie = 0; ie < nEnergies - 1; ie++) {
847 double logE0 = log(generatedE[ireg][ie]);
848 double logE1 = log(generatedE[ireg][ie + 1]);
849 double logBoundary = 0.5 * (logE0 + logE1);
850 boundaryE[ireg][ie] = exp(logBoundary);
851 B2INFO(" " << ie << " " << generatedE[ireg][ie] << " " << boundaryE[ireg][ie] << " GeV");
852 }
853
854 //..Last boundary is an arbitrarily large number
855 boundaryE[ireg][nEnergies - 1] = 9999.;
856 B2INFO(" " << nEnergies - 1 << " " << generatedE[ireg][nEnergies - 1] << " " << boundaryE[ireg][nEnergies - 1] << " GeV");
857 }
858
859 //..Group number of each cellID
860 std::vector<int> groupNumber;
861 for (int cellID = 1; cellID <= ECLElementNumbers::c_NCrystals; cellID++) {
862 const int iGroup = (int)(0.5 + groupNumberOfEachCellID->GetBinContent(cellID));
863 groupNumber.push_back(iGroup);
864 }
865
866 //..Copy results for each group into the payload histogram
867 TH2F nOptimal2D("nOptimal2D", "Optimal nCrys, energy bin vs group number", nCrystalGroups, 0, nCrystalGroups, nEnergies, 0.,
868 nEnergies);
869 for (int ig = 0; ig < nCrystalGroups; ig++) {
870 for (int ie = 0; ie < nEnergies; ie++) {
871 double nOpt = nOptimalPerGroup->GetBinContent(ig + 1, ie + 1);
872 nOptimal2D.SetBinContent(ig + 1, ie + 1, nOpt);
873 }
874 }
875
876 //..2D histograms for peak and bias corrections.
877 // Note that these have 3 x bins per group = nominal energy plus adjacent energies.
878 TH2F peakFracEnergy("peakFracEnergy", "peak contained energy over generated E, energy bin vs group number", 3 * nCrystalGroups, 0,
879 nCrystalGroups, nEnergies, 0., nEnergies);
880 TH2F bias("bias", "median bias (GeV), energy bin vs group number", 3 * nCrystalGroups, 0, nCrystalGroups, nEnergies, 0., nEnergies);
881 TH2F logPeakEnergy("logPeakEnergy", "log peak Energy (GeV), energy bin vs group number", 3 * nCrystalGroups, 0, nCrystalGroups,
882 nEnergies, 0., nEnergies);
883 for (int ig = 0; ig < nCrystalGroups; ig++) {
884 for (int ie = 0; ie < nEnergies; ie++) {
885 const int iy = ie + 1;
886 const int ix = 1 + 3 * ig; // three bins per group in payload histograms
887
888 //..Peak fractional energy and log peak contained energy (GeV)
889
890 //..Nominal ig and and ie
891 double peakAdj = fracContainedEPerGroup->GetBinContent(ig + 1, ie + 1);
892 const double genE = generatedEPerGroup->GetBinContent(ig + 1, ie + 1);
893 double logE = log(genE * peakAdj);
894
895 peakFracEnergy.SetBinContent(ix, iy, peakAdj);
896 logPeakEnergy.SetBinContent(ix, iy, logE);
897
898 //..Lower E adjacent point
899 double peakAdj0 = peakAdj;
900 double logE0 = logE;
901 if (ie > 1) {
902 peakAdj0 = fracContainedAdjacent[0]->GetBinContent(ig + 1, ie + 1);
903 const double genE0 = generatedEPerGroup->GetBinContent(ig + 1, ie);
904 logE0 = log(genE0 * peakAdj0);
905 }
906
907 peakFracEnergy.SetBinContent(ix + 1, iy, peakAdj0);
908 logPeakEnergy.SetBinContent(ix + 1, iy, logE0);
909
910 //..Higher E adjacent point
911 double peakAdj1 = peakAdj;
912 double logE1 = logE;
913 if (ie < nEnergies - 1) {
914 peakAdj1 = fracContainedAdjacent[1]->GetBinContent(ig + 1, ie + 1);
915 const double genE1 = generatedEPerGroup->GetBinContent(ig + 1, ie + 2);
916 logE1 = log(genE1 * peakAdj1);
917 }
918
919 peakFracEnergy.SetBinContent(ix + 2, iy, peakAdj1);
920 logPeakEnergy.SetBinContent(ix + 2, iy, logE1);
921
922 //..Bias
923 double medianBias = biasPerGroup->GetBinContent(ig + 1, ie + 1);
924 bias.SetBinContent(ix, iy, medianBias);
925 medianBias = biasAdjacent[0]->GetBinContent(ig + 1, ie + 1);
926 bias.SetBinContent(ix + 1, iy, medianBias);
927 medianBias = biasAdjacent[1]->GetBinContent(ig + 1, ie + 1);
928 bias.SetBinContent(ix + 2, iy, medianBias);
929 }
930 }
931
932 //..Store the histograms
933 histFile->cd();
934 nOptimal2D.Write();
935 peakFracEnergy.Write();
936 logPeakEnergy.Write();
937 bias.Write();
938 histFile->Close();
939
940 //-----------------------------------------------------------------------------------
941 //..The payload itself
942 ECLnOptimal* optimalNCrys = new ECLnOptimal();
943
944 //..Energy boundaries
945 for (int ireg = 0; ireg < nLeakReg; ireg++) {
946 std::vector<float> eUpperBoundaries;
947 for (int ie = 0; ie < nEnergies; ie++) {eUpperBoundaries.push_back(boundaryE[ireg][ie]);}
948 if (ireg == 0) {
949 optimalNCrys->setUpperBoundariesFwd(eUpperBoundaries);
950 } else if (ireg == 1) {
951 optimalNCrys->setUpperBoundariesBrl(eUpperBoundaries);
952 } else {
953 optimalNCrys->setUpperBoundariesBwd(eUpperBoundaries);
954 }
955 }
956
957 //..Group number of each cellID
958 optimalNCrys->setGroupNumber(groupNumber);
959
960 //..nOptimal histogram
961 optimalNCrys->setNOptimal(nOptimal2D);
962
963 //..Peak and bias correction histograms
964 optimalNCrys->setPeakFracEnergy(peakFracEnergy);
965 optimalNCrys->setBias(bias);
966 optimalNCrys->setLogPeakEnergy(logPeakEnergy);
967
968 //..Save the calibration
969 saveCalibration(optimalNCrys, "ECLnOptimal");
970 B2RESULT("eclNOptimalAlgorithm: successfully stored payload ECLnOptimal");
971
972 //..Done
973 return c_OK;
974}
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
const std::vector< Calibration::ExpRun > & getRunList() const
Get the list of runs for which calibration is called.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Singleton class to cache database objects.
Definition: DBStore.h:31
static DataStore & Instance()
Instance of singleton Store.
Definition: DataStore.cc:54
void setInitializeActive(bool active)
Setter for m_initializeActive.
Definition: DataStore.cc:94
DB object to store the optimal number of crystals to be used in a cluster energy sum,...
Definition: ECLnOptimal.h:23
void setNOptimal(const TH2F &nOptimal)
Set the 2D histogram of nOptimal.
Definition: ECLnOptimal.h:92
void setLogPeakEnergy(const TH2F &logPeakEnergy)
Set the 2D histogram of log of peak contained energy in GeV.
Definition: ECLnOptimal.h:101
void setGroupNumber(const std::vector< int > &groupNumber)
Set the vector of group numbers for each crystal.
Definition: ECLnOptimal.h:89
void setUpperBoundariesFwd(const std::vector< float > &eUpperBoundariesFwd)
Set vector of energies that specify energy ranges in the forward endcap.
Definition: ECLnOptimal.h:80
void setBias(const TH2F &bias)
Set the 2D histogram of beam background bias (reconstructed - true) (GeV)
Definition: ECLnOptimal.h:98
void setUpperBoundariesBrl(const std::vector< float > &eUpperBoundariesBrl)
Set vector of energies that specify energy ranges in the barrel.
Definition: ECLnOptimal.h:83
void setUpperBoundariesBwd(const std::vector< float > &eUpperBoundariesBwd)
Set vector of energies that specify energy ranges in the backward endcap.
Definition: ECLnOptimal.h:86
void setPeakFracEnergy(const TH2F &peakFracEnergy)
Set the 2D histogram of peak fractional contained energy.
Definition: ECLnOptimal.h:95
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
bool isForward(int cellId)
Check whether the crystal is in forward ECL.
const int c_NCrystals
Number of crystals.
bool isBackward(int cellId)
Check whether the crystal is in backward ECL.

◆ 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.

29{
30 // Is it a sequence?
31 if (PySequence_Check(pyObj)) {
32 Py_ssize_t nObj = PySequence_Length(pyObj);
33 // Does it have 2 objects in it?
34 if (nObj != 2) {
35 B2DEBUG(29, "ExpRun was a Python sequence which didn't have exactly 2 entries!");
36 return false;
37 }
38 PyObject* item1, *item2;
39 item1 = PySequence_GetItem(pyObj, 0);
40 item2 = PySequence_GetItem(pyObj, 1);
41 // Did the GetItem work?
42 if ((item1 == NULL) || (item2 == NULL)) {
43 B2DEBUG(29, "A PyObject pointer was NULL in the sequence");
44 return false;
45 }
46 // Are they longs?
47 if (PyLong_Check(item1) && PyLong_Check(item2)) {
48 long value1, value2;
49 value1 = PyLong_AsLong(item1);
50 value2 = PyLong_AsLong(item2);
51 if (((value1 == -1) || (value2 == -1)) && PyErr_Occurred()) {
52 B2DEBUG(29, "An error occurred while converting the PyLong to long");
53 return false;
54 }
55 } else {
56 B2DEBUG(29, "One or more of the PyObjects in the ExpRun wasn't a long");
57 return false;
58 }
59 // Make sure to kill off the reference GetItem gave us responsibility for
60 Py_DECREF(item1);
61 Py_DECREF(item2);
62 } else {
63 B2DEBUG(29, "ExpRun was not a Python sequence.");
64 return false;
65 }
66 return true;
67}

◆ clearCalibrationData()

void clearCalibrationData ( )
inlineprotectedinherited

Clear calibration data.

Definition at line 324 of file CalibrationAlgorithm.h.

void clearCalibrationData()
Clear calibration data.
ExecutionData m_data
Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.

◆ commit() [1/2]

bool commit ( )
inherited

Submit constants from last calibration into database.

Definition at line 302 of file CalibrationAlgorithm.cc.

303{
304 if (getPayloads().empty())
305 return false;
306 list<Database::DBImportQuery> payloads = getPayloads();
307 B2INFO("Committing " << payloads.size() << " payloads to database.");
308 return Database::Instance().storeData(payloads);
309}
std::list< Database::DBImportQuery > & getPayloads()
Get constants (in TObjects) for database update from last execution.
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:42
bool storeData(const std::string &name, TObject *object, const IntervalOfValidity &iov)
Store an object in the database.
Definition: Database.cc:141

◆ commit() [2/2]

bool commit ( std::list< Database::DBImportQuery payloads)
inherited

Submit constants from a (potentially previous) set of payloads.

Definition at line 311 of file CalibrationAlgorithm.cc.

312{
313 if (payloads.empty())
314 return false;
315 return Database::Instance().storeData(payloads);
316}

◆ 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.

71{
72 ExpRun expRun;
73 PyObject* itemExp, *itemRun;
74 itemExp = PySequence_GetItem(pyObj, 0);
75 itemRun = PySequence_GetItem(pyObj, 1);
76 expRun.first = PyLong_AsLong(itemExp);
77 Py_DECREF(itemExp);
78 expRun.second = PyLong_AsLong(itemRun);
79 Py_DECREF(itemRun);
80 return expRun;
81}
Struct containing exp number and run number.
Definition: Splitter.h:51

◆ dumpOutputJson()

const std::string dumpOutputJson ( ) const
inlineinherited

Dump the JSON string of the output JSON object.

Definition at line 223 of file CalibrationAlgorithm.h.

223{return m_jsonExecutionOutput.dump();}
nlohmann::json m_jsonExecutionOutput
Optional output JSON object that can be set during the execution by the underlying algorithm code.

◆ execute() [1/2]

CalibrationAlgorithm::EResult execute ( PyObject *  runs,
int  iteration = 0,
IntervalOfValidity  iov = IntervalOfValidity() 
)
inherited

Runs calibration over Python list of runs. Converts to C++ and then calls the other execute() function.

Definition at line 83 of file CalibrationAlgorithm.cc.

84{
85 B2DEBUG(29, "Running execute() using Python Object as input argument");
86 // Reset the execution specific data in case the algorithm was previously called
87 m_data.reset();
88 m_data.setIteration(iteration);
89 vector<ExpRun> vecRuns;
90 // Is it a list?
91 if (PySequence_Check(runs)) {
92 boost::python::handle<> handle(boost::python::borrowed(runs));
93 boost::python::list listRuns(handle);
94
95 int nList = boost::python::len(listRuns);
96 for (int iList = 0; iList < nList; ++iList) {
97 boost::python::object pyExpRun(listRuns[iList]);
98 if (!checkPyExpRun(pyExpRun.ptr())) {
99 B2ERROR("Received Python ExpRuns couldn't be converted to C++");
101 return c_Failure;
102 } else {
103 vecRuns.push_back(convertPyExpRun(pyExpRun.ptr()));
104 }
105 }
106 } else {
107 B2ERROR("Tried to set the input runs but we didn't receive a Python sequence object (list,tuple).");
109 return c_Failure;
110 }
111 return execute(vecRuns, iteration, iov);
112}
void setResult(EResult result)
Setter for current iteration.
void setIteration(int iteration)
Setter for current iteration.
void reset()
Resets this class back to what is needed at the beginning of an execution.
bool checkPyExpRun(PyObject *pyObj)
Checks that a PyObject can be successfully converted to an ExpRun type.
EResult execute(std::vector< Calibration::ExpRun > runs={}, int iteration=0, IntervalOfValidity iov=IntervalOfValidity())
Runs calibration over vector of runs for a given iteration.
Calibration::ExpRun convertPyExpRun(PyObject *pyObj)
Performs the conversion of PyObject to ExpRun.

◆ execute() [2/2]

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.

115{
116 // Check if we are calling this function directly and need to reset, or through Python where it was already done.
117 if (m_data.getResult() != c_Undefined) {
118 m_data.reset();
119 m_data.setIteration(iteration);
120 }
121
122 if (m_inputFileNames.empty()) {
123 B2ERROR("There aren't any input files set. Please use CalibrationAlgorithm::setInputFiles()");
125 return c_Failure;
126 }
127
128 // Did we receive runs to execute over explicitly?
129 if (!(runs.empty())) {
130 for (auto expRun : runs) {
131 B2DEBUG(29, "ExpRun requested = (" << expRun.first << ", " << expRun.second << ")");
132 }
133 // We've asked explicitly for certain runs, but we should check if the data granularity is 'run'
134 if (strcmp(getGranularity().c_str(), "all") == 0) {
135 B2ERROR(("The data is collected with granularity=all (exp=-1,run=-1), but you seem to request calibration for specific runs."
136 " We'll continue but using ALL the input data given instead of the specific runs requested."));
137 }
138 } else {
139 // If no runs are provided, infer the runs from all collected data
140 runs = getRunListFromAllData();
141 // Let's check that we have some now
142 if (runs.empty()) {
143 B2ERROR("No collected data in input files.");
145 return c_Failure;
146 }
147 for (auto expRun : runs) {
148 B2DEBUG(29, "ExpRun requested = (" << expRun.first << ", " << expRun.second << ")");
149 }
150 }
151
153 if (iov.empty()) {
154 // If no user specified IoV we use the IoV from the executed run list
155 iov = IntervalOfValidity(runs[0].first, runs[0].second, runs[runs.size() - 1].first, runs[runs.size() - 1].second);
156 }
158 // After here, the getObject<...>(...) helpers start to work
159
161 m_data.setResult(result);
162 return result;
163}
void setRequestedIov(const IntervalOfValidity &iov=IntervalOfValidity(0, 0, -1, -1))
Sets the requested IoV for this execution, based on the.
void setRequestedRuns(const std::vector< Calibration::ExpRun > &requestedRuns)
Sets the vector of ExpRuns.
EResult getResult() const
Getter for current result.
std::vector< Calibration::ExpRun > getRunListFromAllData() const
Get the complete list of runs from inspection of collected data.
std::vector< std::string > m_inputFileNames
List of input files to the Algorithm, will initially be user defined but then gets the wildcards expa...
EResult
The result of calibration.
virtual EResult calibrate()=0
Run algo on data - pure virtual: needs to be implemented.
std::string getGranularity() const
Get the granularity of collected data.
A class that describes the interval of experiments/runs for which an object in the database is valid.

◆ fillRunToInputFilesMap()

void fillRunToInputFilesMap ( )
inherited

Fill the mapping of ExpRun -> Files.

Definition at line 330 of file CalibrationAlgorithm.cc.

331{
332 m_runsToInputFiles.clear();
333 // Save TDirectory to change back at the end
334 TDirectory* dir = gDirectory;
335 RunRange* runRange;
336 // Construct the TDirectory name where we expect our objects to be
337 string runRangeObjName(getPrefix() + "/" + RUN_RANGE_OBJ_NAME);
338 for (const auto& fileName : m_inputFileNames) {
339 //Open TFile to get the objects
340 unique_ptr<TFile> f;
341 f.reset(TFile::Open(fileName.c_str(), "READ"));
342 runRange = dynamic_cast<RunRange*>(f->Get(runRangeObjName.c_str()));
343 if (runRange) {
344 // Insert or extend the run -> file mapping for this ExpRun
345 auto expRuns = runRange->getExpRunSet();
346 for (const auto& expRun : expRuns) {
347 auto runFiles = m_runsToInputFiles.find(expRun);
348 if (runFiles != m_runsToInputFiles.end()) {
349 (runFiles->second).push_back(fileName);
350 } else {
351 m_runsToInputFiles.insert(std::make_pair(expRun, std::vector<std::string> {fileName}));
352 }
353 }
354 } else {
355 B2WARNING("Missing a RunRange object for file: " << fileName);
356 }
357 }
358 dir->cd();
359}
std::string getPrefix() const
Get the prefix used for getting calibration data.
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 setti...
Mergeable object holding (unique) set of (exp,run) pairs.
Definition: RunRange.h:25
const std::set< Calibration::ExpRun > & getExpRunSet()
Get access to the stored set.
Definition: RunRange.h:64

◆ findPayloadBoundaries()

const std::vector< ExpRun > findPayloadBoundaries ( std::vector< Calibration::ExpRun >  runs,
int  iteration = 0 
)
inherited

Used to discover the ExpRun boundaries that you want the Python CAF to execute on. This is optional and only used in some.

Definition at line 520 of file CalibrationAlgorithm.cc.

521{
522 m_boundaries.clear();
523 if (m_inputFileNames.empty()) {
524 B2ERROR("There aren't any input files set. Please use CalibrationAlgorithm::setInputFiles()");
525 return m_boundaries;
526 }
527 // Reset the internal execution data just in case something is hanging around
528 m_data.reset();
529 if (runs.empty()) {
530 // Want to loop over all runs we could possibly know about
531 runs = getRunListFromAllData();
532 }
533 // Let's check that we have some now
534 if (runs.empty()) {
535 B2ERROR("No collected data in input files.");
536 return m_boundaries;
537 }
538 // In order to find run boundaries we must have collected with data granularity == 'run'
539 if (strcmp(getGranularity().c_str(), "all") == 0) {
540 B2ERROR("The data is collected with granularity='all' (exp=-1,run=-1), and we can't use that to find run boundaries.");
541 return m_boundaries;
542 }
543 m_data.setIteration(iteration);
544 // User defined setup function
545 boundaryFindingSetup(runs, iteration);
546 std::vector<ExpRun> runList;
547 // Loop over run list and call derived class "isBoundaryRequired" member function
548 for (auto currentRun : runs) {
549 runList.push_back(currentRun);
550 m_data.setRequestedRuns(runList);
551 // After here, the getObject<...>(...) helpers start to work
552 if (isBoundaryRequired(currentRun)) {
553 m_boundaries.push_back(currentRun);
554 }
555 // Only want run-by-run
556 runList.clear();
557 // Don't want memory hanging around
559 }
560 m_data.reset();
562 return m_boundaries;
563}
std::vector< Calibration::ExpRun > m_boundaries
When using the boundaries functionality from isBoundaryRequired, this is used to store the boundaries...
virtual void boundaryFindingTearDown()
Put your algorithm back into a state ready for normal execution if you need to.
virtual void boundaryFindingSetup(std::vector< Calibration::ExpRun >, int)
If you need to make some changes to your algorithm class before 'findPayloadBoundaries' is run,...
virtual bool isBoundaryRequired(const Calibration::ExpRun &)
Given the current collector data, make a decision about whether or not this run should be the start o...

◆ getAllGranularityExpRun()

Calibration::ExpRun getAllGranularityExpRun ( ) const
inlineprotectedinherited

Returns the Exp,Run pair that means 'Everything'. Currently unused.

Definition at line 327 of file CalibrationAlgorithm.h.

327{return m_allExpRun;}
static const Calibration::ExpRun m_allExpRun
allExpRun

◆ 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 algorithm.

Definition at line 164 of file CalibrationAlgorithm.h.

164{return getPrefix();}

◆ getDescription()

const std::string & getDescription ( ) const
inlineinherited

Get the description of the algorithm (set by developers in constructor)

Definition at line 216 of file CalibrationAlgorithm.h.

216{return m_description;}
std::string m_description
Description of the algorithm.

◆ getExpRunString()

string getExpRunString ( Calibration::ExpRun &  expRun) const
privateinherited

Gets the "exp.run" string repr. of (exp,run)

Definition at line 254 of file CalibrationAlgorithm.cc.

255{
256 string expRunString;
257 expRunString += to_string(expRun.first);
258 expRunString += ".";
259 expRunString += to_string(expRun.second);
260 return expRunString;
261}

◆ getFullObjectPath()

string getFullObjectPath ( const std::string &  name,
Calibration::ExpRun  expRun 
) const
privateinherited

constructs the full TDirectory + Key name of an object in a TFile based on its name and exprun

Definition at line 263 of file CalibrationAlgorithm.cc.

264{
265 string dirName = getPrefix() + "/" + name;
266 string objName = name + "_" + getExpRunString(expRun);
267 return dirName + "/" + objName;
268}
std::string getExpRunString(Calibration::ExpRun &expRun) const
Gets the "exp.run" string repr. of (exp,run)

◆ getGranularity()

std::string getGranularity ( ) const
inlineinherited

Get the granularity of collected data.

Definition at line 188 of file CalibrationAlgorithm.h.

188{return m_granularityOfData;};
std::string m_granularityOfData
Granularity of input data. This only changes when the input files change so it isn't specific to an e...

◆ getGranularityFromData()

string getGranularityFromData ( ) const
protectedinherited

Get the granularity of collected data.

Definition at line 383 of file CalibrationAlgorithm.cc.

384{
385 // Save TDirectory to change back at the end
386 TDirectory* dir = gDirectory;
387 RunRange* runRange;
388 string runRangeObjName(getPrefix() + "/" + RUN_RANGE_OBJ_NAME);
389 // We only check the first file
390 string fileName = m_inputFileNames[0];
391 unique_ptr<TFile> f;
392 f.reset(TFile::Open(fileName.c_str(), "READ"));
393 runRange = dynamic_cast<RunRange*>(f->Get(runRangeObjName.c_str()));
394 if (!runRange) {
395 B2FATAL("The input file " << fileName << " does not contain a RunRange object at "
396 << runRangeObjName << ". Please set your input files to exclude it.");
397 return "";
398 }
399 string granularity = runRange->getGranularity();
400 dir->cd();
401 return granularity;
402}
std::string getGranularity() const
Gets the m_granularity.
Definition: RunRange.h:110

◆ getInputFileNames()

PyObject * getInputFileNames ( )
inherited

Get the input file names used for this algorithm and pass them out as a Python list of unicode strings.

Definition at line 245 of file CalibrationAlgorithm.cc.

246{
247 PyObject* objInputFileNames = PyList_New(m_inputFileNames.size());
248 for (size_t i = 0; i < m_inputFileNames.size(); ++i) {
249 PyList_SetItem(objInputFileNames, i, Py_BuildValue("s", m_inputFileNames[i].c_str()));
250 }
251 return objInputFileNames;
252}

◆ getInputJsonObject()

const nlohmann::json & getInputJsonObject ( ) const
inlineprotectedinherited

Get the entire top level JSON object. We explicitly say this must be of object type so that we might pick.

Definition at line 357 of file CalibrationAlgorithm.h.

357{return m_jsonExecutionInput;}
nlohmann::json m_jsonExecutionInput
Optional input JSON object used to make decisions about how to execute the algorithm code.

◆ getInputJsonValue()

const T getInputJsonValue ( const std::string &  key) const
inlineprotectedinherited

Get an input JSON value using a key. The normal exceptions are raised when the key doesn't exist.

Definition at line 350 of file CalibrationAlgorithm.h.

351 {
352 return m_jsonExecutionInput.at(key);
353 }

◆ getIovFromAllData()

IntervalOfValidity getIovFromAllData ( ) const
inherited

Get the complete IoV from inspection of collected data.

Definition at line 325 of file CalibrationAlgorithm.cc.

326{
328}
RunRange getRunRangeFromAllData() const
Get the complete RunRange from inspection of collected data.
IntervalOfValidity getIntervalOfValidity()
Make IntervalOfValidity from the set, spanning all runs. Works because sets are sorted by default.
Definition: RunRange.h:70

◆ getIteration()

int getIteration ( ) const
inlineprotectedinherited

Get current iteration.

Definition at line 269 of file CalibrationAlgorithm.h.

269{ return m_data.getIteration(); }
int getIteration() const
Getter for current iteration.

◆ getObjectPtr()

std::shared_ptr< T > getObjectPtr ( std::string  name)
inlineprotectedinherited

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.

Definition at line 285 of file CalibrationAlgorithm.h.

286 {
287 if (m_runsToInputFiles.size() == 0)
289 return getObjectPtr<T>(name, m_data.getRequestedRuns());
290 }
const std::vector< Calibration::ExpRun > & getRequestedRuns() const
Returns the vector of ExpRuns.
void fillRunToInputFilesMap()
Fill the mapping of ExpRun -> Files.

◆ getOutputJsonValue()

const T getOutputJsonValue ( const std::string &  key) const
inlineprotectedinherited

Get a value using a key from the JSON output object, not sure why you would want to do this.

Definition at line 342 of file CalibrationAlgorithm.h.

343 {
344 return m_jsonExecutionOutput.at(key);
345 }

◆ getPayloads()

std::list< Database::DBImportQuery > & getPayloads ( )
inlineinherited

Get constants (in TObjects) for database update from last execution.

Definition at line 204 of file CalibrationAlgorithm.h.

204{return m_data.getPayloads();}
std::list< Database::DBImportQuery > & getPayloads()
Get constants (in TObjects) for database update from last calibration.

◆ getPayloadValues()

std::list< Database::DBImportQuery > getPayloadValues ( )
inlineinherited

Get constants (in TObjects) for database update from last execution but passed by VALUE.

Definition at line 207 of file CalibrationAlgorithm.h.

207{return m_data.getPayloadValues();}
std::list< Database::DBImportQuery > getPayloadValues()
Get constants (in TObjects) for database update from last calibration but passed by VALUE.

◆ getPrefix()

std::string getPrefix ( ) const
inlineinherited

Get the prefix used for getting calibration data.

Definition at line 146 of file CalibrationAlgorithm.h.

146{return m_prefix;}
std::string m_prefix
The name of the TDirectory the collector objects are contained within.

◆ getRunList()

const std::vector< Calibration::ExpRun > & getRunList ( ) const
inlineprotectedinherited

Get the list of runs for which calibration is called.

Definition at line 266 of file CalibrationAlgorithm.h.

266{return m_data.getRequestedRuns();}

◆ getRunListFromAllData()

vector< ExpRun > getRunListFromAllData ( ) const
inherited

Get the complete list of runs from inspection of collected data.

Definition at line 318 of file CalibrationAlgorithm.cc.

319{
320 RunRange runRange = getRunRangeFromAllData();
321 set<ExpRun> expRunSet = runRange.getExpRunSet();
322 return vector<ExpRun>(expRunSet.begin(), expRunSet.end());
323}

◆ getRunRangeFromAllData()

RunRange getRunRangeFromAllData ( ) const
inherited

Get the complete RunRange from inspection of collected data.

Definition at line 361 of file CalibrationAlgorithm.cc.

362{
363 // Save TDirectory to change back at the end
364 TDirectory* dir = gDirectory;
365 RunRange runRange;
366 // Construct the TDirectory name where we expect our objects to be
367 string runRangeObjName(getPrefix() + "/" + RUN_RANGE_OBJ_NAME);
368 for (const auto& fileName : m_inputFileNames) {
369 //Open TFile to get the objects
370 unique_ptr<TFile> f;
371 f.reset(TFile::Open(fileName.c_str(), "READ"));
372 RunRange* runRangeOther = dynamic_cast<RunRange*>(f->Get(runRangeObjName.c_str()));
373 if (runRangeOther) {
374 runRange.merge(runRangeOther);
375 } else {
376 B2WARNING("Missing a RunRange object for file: " << fileName);
377 }
378 }
379 dir->cd();
380 return runRange;
381}
virtual void merge(const RunRange *other)
Implementation of merging - other is added to the set (union)
Definition: RunRange.h:52

◆ getVecInputFileNames()

std::vector< std::string > getVecInputFileNames ( ) const
inlineprotectedinherited

Get the input file names used for this algorithm as a STL vector.

Definition at line 275 of file CalibrationAlgorithm.h.

275{return m_inputFileNames;}

◆ inputJsonKeyExists()

bool inputJsonKeyExists ( const std::string &  key) const
inlineprotectedinherited

Test for a key in the input JSON object.

Definition at line 360 of file CalibrationAlgorithm.h.

360{return m_jsonExecutionInput.count(key);}

◆ isBoundaryRequired()

virtual bool isBoundaryRequired ( const Calibration::ExpRun &  )
inlineprotectedvirtualinherited

Given the current collector data, make a decision about whether or not this run should be the start of a payload boundary.

Reimplemented in TestBoundarySettingAlgorithm, PXDAnalyticGainCalibrationAlgorithm, PXDValidationAlgorithm, TestCalibrationAlgorithm, SVD3SampleCoGTimeCalibrationAlgorithm, SVD3SampleELSTimeCalibrationAlgorithm, and SVDCoGTimeCalibrationAlgorithm.

Definition at line 243 of file CalibrationAlgorithm.h.

244 {
245 B2ERROR("You didn't implement a isBoundaryRequired() member function in your CalibrationAlgorithm but you are calling it!");
246 return false;
247 }

◆ loadInputJson()

bool loadInputJson ( const std::string &  jsonString)
inherited

Load the m_inputJson variable from a string (useful from Python interface). The return bool indicates success or failure.

Definition at line 502 of file CalibrationAlgorithm.cc.

503{
504 try {
505 auto jsonInput = nlohmann::json::parse(jsonString);
506 // Input string has an object (dict) as the top level object?
507 if (jsonInput.is_object()) {
508 m_jsonExecutionInput = jsonInput;
509 return true;
510 } else {
511 B2ERROR("JSON input string isn't an object type i.e. not a '{}' at the top level.");
512 return false;
513 }
514 } catch (nlohmann::json::parse_error&) {
515 B2ERROR("Parsing of JSON input string failed");
516 return false;
517 }
518}

◆ resetInputJson()

void resetInputJson ( )
inlineprotectedinherited

Clears the m_inputJson member variable.

Definition at line 330 of file CalibrationAlgorithm.h.

330{m_jsonExecutionInput.clear();}

◆ resetOutputJson()

void resetOutputJson ( )
inlineprotectedinherited

Clears the m_outputJson member variable.

Definition at line 333 of file CalibrationAlgorithm.h.

333{m_jsonExecutionOutput.clear();}

◆ saveCalibration() [1/6]

void saveCalibration ( TClonesArray *  data,
const std::string &  name 
)
protectedinherited

Store DBArray payload with given name with default IOV.

Definition at line 297 of file CalibrationAlgorithm.cc.

298{
300}
const IntervalOfValidity & getRequestedIov() const
Getter for requested IOV.

◆ saveCalibration() [2/6]

void saveCalibration ( TClonesArray *  data,
const std::string &  name,
const IntervalOfValidity iov 
)
protectedinherited

Store DBArray with given name and custom IOV.

Definition at line 276 of file CalibrationAlgorithm.cc.

277{
278 B2DEBUG(29, "Saving calibration TClonesArray '" << name << "' to payloads list.");
279 getPayloads().emplace_back(name, data, iov);
280}

◆ saveCalibration() [3/6]

void saveCalibration ( TObject *  data)
protectedinherited

Store DB payload with default name and default IOV.

Definition at line 287 of file CalibrationAlgorithm.cc.

288{
289 saveCalibration(data, DataStore::objectName(data->IsA(), ""));
290}
static std::string objectName(const TClass *t, const std::string &name)
Return the storage name for an object of the given TClass and name.
Definition: DataStore.cc:151

◆ saveCalibration() [4/6]

void saveCalibration ( TObject *  data,
const IntervalOfValidity iov 
)
protectedinherited

Store DB payload with default name and custom IOV.

Definition at line 282 of file CalibrationAlgorithm.cc.

283{
284 saveCalibration(data, DataStore::objectName(data->IsA(), ""), iov);
285}

◆ saveCalibration() [5/6]

void saveCalibration ( TObject *  data,
const std::string &  name 
)
protectedinherited

Store DB payload with given name with default IOV.

Definition at line 292 of file CalibrationAlgorithm.cc.

293{
295}

◆ saveCalibration() [6/6]

void saveCalibration ( TObject *  data,
const std::string &  name,
const IntervalOfValidity iov 
)
protectedinherited

Store DB payload with given name and custom IOV.

Definition at line 270 of file CalibrationAlgorithm.cc.

271{
272 B2DEBUG(29, "Saving calibration TObject = '" << name << "' to payloads list.");
273 getPayloads().emplace_back(name, data, iov);
274}

◆ setDescription()

void setDescription ( const std::string &  description)
inlineprotectedinherited

Set algorithm description (in constructor)

Definition at line 321 of file CalibrationAlgorithm.h.

321{m_description = description;}

◆ 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.

167{
168 // The reasoning for this very 'manual' approach to extending the Python interface
169 // (instead of using boost::python) is down to my fear of putting off final users with
170 // complexity on their side.
171 //
172 // I didn't want users that inherit from this class to be forced to use boost and
173 // to have to define a new python module just to use the CAF. A derived class from
174 // from a boost exposed class would need to have its own boost python module definition
175 // to allow access from a steering file and to the base class functions (I think).
176 // I also couldn't be bothered to write a full framework to get around the issue in a similar
177 // way to Module()...maybe there's an easy way.
178 //
179 // But this way we can allow people to continue using their ROOT implemented classes and inherit
180 // easily from this one. But add in a few helper functions that work with Python objects
181 // created in their steering file i.e. instead of being forced to use STL objects as input
182 // to the algorithm.
183 if (PyList_Check(inputFileNames)) {
184 boost::python::handle<> handle(boost::python::borrowed(inputFileNames));
185 boost::python::list listInputFileNames(handle);
186 auto vecInputFileNames = PyObjConvUtils::convertPythonObject(listInputFileNames, vector<string>());
187 setInputFileNames(vecInputFileNames);
188 } else {
189 B2ERROR("Tried to set the input files but we didn't receive a Python list.");
190 }
191}
void setInputFileNames(PyObject *inputFileNames)
Set the input file names used for this algorithm from a Python list.
Scalar convertPythonObject(const boost::python::object &pyObject, Scalar)
Convert from Python to given type.

◆ 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.

195{
196 // A lot of code below is tweaked from RootInputModule::initialize,
197 // since we're basically copying the functionality anyway.
198 if (inputFileNames.empty()) {
199 B2WARNING("You have called setInputFileNames() with an empty list. Did you mean to do that?");
200 return;
201 }
202 auto tmpInputFileNames = RootIOUtilities::expandWordExpansions(inputFileNames);
203
204 // We'll use a set to enforce sorted unique file paths as we check them
205 set<string> setInputFileNames;
206 // Check that files exist and convert to absolute paths
207 for (auto path : tmpInputFileNames) {
208 string fullPath = fs::absolute(path).string();
209 if (fs::exists(fullPath)) {
210 setInputFileNames.insert(fs::canonical(fullPath).string());
211 } else {
212 B2WARNING("Couldn't find the file " << path);
213 }
214 }
215
216 if (setInputFileNames.empty()) {
217 B2WARNING("No valid files specified!");
218 return;
219 } else {
220 // Reset the run -> files map as our files are likely different
221 m_runsToInputFiles.clear();
222 }
223
224 // Open TFile to check they can be accessed by ROOT
225 TDirectory* dir = gDirectory;
226 for (const string& fileName : setInputFileNames) {
227 unique_ptr<TFile> f;
228 try {
229 f.reset(TFile::Open(fileName.c_str(), "READ"));
230 } catch (logic_error&) {
231 //this might happen for ~invaliduser/foo.root
232 //actually undefined behaviour per standard, reported as ROOT-8490 in JIRA
233 }
234 if (!f || !f->IsOpen()) {
235 B2FATAL("Couldn't open input file " + fileName);
236 }
237 }
238 dir->cd();
239
240 // Copy the entries of the set to a vector
241 m_inputFileNames = vector<string>(setInputFileNames.begin(), setInputFileNames.end());
243}
std::string getGranularityFromData() const
Get the granularity of collected data.
std::vector< std::string > expandWordExpansions(const std::vector< std::string > &filenames)
Performs wildcard expansion using wordexp(), returns matches.

◆ setOutputJsonValue()

void setOutputJsonValue ( const std::string &  key,
const T &  value 
)
inlineprotectedinherited

Set a key:value pair for the outputJson object, expected to used internally during calibrate()

Definition at line 337 of file CalibrationAlgorithm.h.

337{m_jsonExecutionOutput[key] = value;}

◆ setPrefix()

void setPrefix ( const std::string &  prefix)
inlineinherited

Set the prefix used to identify datastore objects.

Definition at line 167 of file CalibrationAlgorithm.h.

167{m_prefix = prefix;}

◆ updateDBObjPtrs()

void updateDBObjPtrs ( const unsigned int  event = 1,
const int  run = 0,
const int  experiment = 0 
)
protectedinherited

Updates any DBObjPtrs by calling update(event) for DBStore.

Definition at line 404 of file CalibrationAlgorithm.cc.

405{
406 // Construct an EventMetaData object but NOT in the Datastore
407 EventMetaData emd(event, run, experiment);
408 // Explicitly update while avoiding registering a Datastore object
410 // Also update the intra-run objects to the event at the same time (maybe unnecessary...)
412}
Store event, run, and experiment numbers.
Definition: EventMetaData.h:33

Member Data Documentation

◆ m_allExpRun

const ExpRun m_allExpRun = make_pair(-1, -1)
staticprivateinherited

allExpRun

Definition at line 364 of file CalibrationAlgorithm.h.

◆ m_boundaries

std::vector<Calibration::ExpRun> m_boundaries
protectedinherited

When using the boundaries functionality from isBoundaryRequired, this is used to store the boundaries. It is cleared when.

Definition at line 261 of file CalibrationAlgorithm.h.

◆ m_data

ExecutionData m_data
privateinherited

Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.

Definition at line 382 of file CalibrationAlgorithm.h.

◆ m_description

std::string m_description {""}
privateinherited

Description of the algorithm.

Definition at line 385 of file CalibrationAlgorithm.h.

◆ m_granularityOfData

std::string m_granularityOfData
privateinherited

Granularity of input data. This only changes when the input files change so it isn't specific to an execution.

Definition at line 379 of file CalibrationAlgorithm.h.

◆ m_inputFileNames

std::vector<std::string> m_inputFileNames
privateinherited

List of input files to the Algorithm, will initially be user defined but then gets the wildcards expanded during execute()

Definition at line 373 of file CalibrationAlgorithm.h.

◆ m_jsonExecutionInput

nlohmann::json m_jsonExecutionInput = nlohmann::json::object()
privateinherited

Optional input JSON object used to make decisions about how to execute the algorithm code.

Definition at line 397 of file CalibrationAlgorithm.h.

◆ m_jsonExecutionOutput

nlohmann::json m_jsonExecutionOutput = nlohmann::json::object()
privateinherited

Optional output JSON object that can be set during the execution by the underlying algorithm code.

Definition at line 403 of file CalibrationAlgorithm.h.

◆ m_prefix

std::string m_prefix {""}
privateinherited

The name of the TDirectory the collector objects are contained within.

Definition at line 388 of file CalibrationAlgorithm.h.

◆ m_runsToInputFiles

std::map<Calibration::ExpRun, std::vector<std::string> > m_runsToInputFiles
privateinherited

Map of Runs to input files. Gets filled when you call getRunRangeFromAllData, gets cleared when setting input files again.

Definition at line 376 of file CalibrationAlgorithm.h.


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