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
 
const std::string & getPrefix () const
 Get the prefix used for getting calibration data.
 
const 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.
 
const 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.
 
const std::list< Database::DBImportQuery > & getPayloadValues () const
 Get constants (in TObjects) for database update from last execution.
 
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>.
 

Static Public Member Functions

static bool checkPyExpRun (PyObject *pyObj)
 Checks that a PyObject can be successfully converted to an ExpRun type.
 
static Calibration::ExpRun convertPyExpRun (PyObject *pyObj)
 Performs the conversion of PyObject to ExpRun.
 

Protected Member Functions

virtual EResult calibrate () override
 ..Run algorithm on events
 
void setInputFileNames (const 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.
 
const 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 setDescription (const std::string &description)
 Set algorithm description (in constructor)
 
void clearCalibrationData ()
 Clear calibration data.
 
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.
 

Static Protected Member Functions

static void updateDBObjPtrs (const unsigned int event, const int run, const int experiment)
 Updates any DBObjPtrs by calling update(event) for DBStore.
 
static Calibration::ExpRun getAllGranularityExpRun ()
 Returns the Exp,Run pair that means 'Everything'. Currently unused.
 

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,
43 c_NotEnoughData,
44 c_Failure,
45 c_Undefined
46 };

Constructor & Destructor Documentation

◆ eclNOptimalAlgorithm()

..Constructor

Definition at line 223 of file eclNOptimalAlgorithm.cc.

223 : CalibrationAlgorithm("eclNOptimalCollector")
224{
226 "Use single gamma MC to find the optimal number of crystals to sum for cluster energy"
227 );
228}
void setDescription(const std::string &description)
Set algorithm description (in constructor)
CalibrationAlgorithm(const std::string &collectorModuleName)
Constructor - sets the prefix for collected objects (won't be accesses until execute(....

◆ ~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 PXDAnalyticGainCalibrationAlgorithm, PXDValidationAlgorithm, SVD3SampleCoGTimeCalibrationAlgorithm, SVD3SampleELSTimeCalibrationAlgorithm, SVDClusterAbsoluteTimeShifterAlgorithm, SVDCoGTimeCalibrationAlgorithm, TestBoundarySettingAlgorithm, and TestCalibrationAlgorithm.

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 232 of file eclNOptimalAlgorithm.cc.

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

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.

324{m_data.clearCalibrationData();}

◆ 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 312 of file CalibrationAlgorithm.cc.

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

◆ convertPyExpRun()

ExpRun convertPyExpRun ( PyObject * pyObj)
staticinherited

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}

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

◆ 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++");
100 m_data.setResult(c_Failure);
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).");
108 m_data.setResult(c_Failure);
109 return c_Failure;
110 }
111 return execute(vecRuns, iteration, iov);
112}
static 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.
static Calibration::ExpRun convertPyExpRun(PyObject *pyObj)
Performs the conversion of PyObject to ExpRun.
ExecutionData m_data
Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.

◆ 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()");
124 m_data.setResult(c_Failure);
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.");
144 m_data.setResult(c_Failure);
145 return c_Failure;
146 }
147 for (auto expRun : runs) {
148 B2DEBUG(29, "ExpRun requested = (" << expRun.first << ", " << expRun.second << ")");
149 }
150 }
151
152 m_data.setRequestedRuns(runs);
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 }
157 m_data.setRequestedIov(iov);
158 // After here, the getObject<...>(...) helpers start to work
159
161 m_data.setResult(result);
162 return result;
163}
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.
@ c_Undefined
Not yet known (before execution) =4 in Python.
const std::string & getGranularity() const
Get the granularity of collected data.
virtual EResult calibrate()=0
Run algo on data - pure virtual: needs to be implemented.

◆ fillRunToInputFilesMap()

void fillRunToInputFilesMap ( )
inherited

Fill the mapping of ExpRun -> Files.

Definition at line 331 of file CalibrationAlgorithm.cc.

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

522{
523 m_boundaries.clear();
524 if (m_inputFileNames.empty()) {
525 B2ERROR("There aren't any input files set. Please use CalibrationAlgorithm::setInputFiles()");
526 return m_boundaries;
527 }
528 // Reset the internal execution data just in case something is hanging around
529 m_data.reset();
530 if (runs.empty()) {
531 // Want to loop over all runs we could possibly know about
532 runs = getRunListFromAllData();
533 }
534 // Let's check that we have some now
535 if (runs.empty()) {
536 B2ERROR("No collected data in input files.");
537 return m_boundaries;
538 }
539 // In order to find run boundaries we must have collected with data granularity == 'run'
540 if (strcmp(getGranularity().c_str(), "all") == 0) {
541 B2ERROR("The data is collected with granularity='all' (exp=-1,run=-1), and we can't use that to find run boundaries.");
542 return m_boundaries;
543 }
544 m_data.setIteration(iteration);
545 // User defined setup function
546 boundaryFindingSetup(runs, iteration);
547 std::vector<ExpRun> runList;
548 // Loop over run list and call derived class "isBoundaryRequired" member function
549 for (auto currentRun : runs) {
550 runList.push_back(currentRun);
551 m_data.setRequestedRuns(runList);
552 // After here, the getObject<...>(...) helpers start to work
553 if (isBoundaryRequired(currentRun)) {
554 m_boundaries.push_back(currentRun);
555 }
556 // Only want run-by-run
557 runList.clear();
558 // Don't want memory hanging around
559 m_data.clearCalibrationData();
560 }
561 m_data.reset();
563 return m_boundaries;
564}
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()

static Calibration::ExpRun getAllGranularityExpRun ( )
inlinestaticprotectedinherited

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

Definition at line 327 of file CalibrationAlgorithm.h.

327{return m_allExpRun;}

◆ getCollectorName()

const 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;}

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

const std::string & getGranularity ( ) const
inlineinherited

Get the granularity of collected data.

Definition at line 188 of file CalibrationAlgorithm.h.

188{return m_granularityOfData;};

◆ getGranularityFromData()

string getGranularityFromData ( ) const
protectedinherited

Get the granularity of collected data.

Definition at line 384 of file CalibrationAlgorithm.cc.

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

◆ getInputJsonValue()

template<class T>
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 326 of file CalibrationAlgorithm.cc.

327{
329}
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(); }

◆ getObjectPtr()

template<class T>
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)
288 fillRunToInputFilesMap();
289 return getObjectPtr<T>(name, m_data.getRequestedRuns());
290 }

◆ getOutputJsonValue()

template<class T>
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();}

◆ getPayloadValues()

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

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

Definition at line 207 of file CalibrationAlgorithm.h.

207{return m_data.getPayloadValues();}

◆ getPrefix()

const 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;}

◆ 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 319 of file CalibrationAlgorithm.cc.

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

◆ getRunRangeFromAllData()

RunRange getRunRangeFromAllData ( ) const
inherited

Get the complete RunRange from inspection of collected data.

Definition at line 362 of file CalibrationAlgorithm.cc.

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

◆ getVecInputFileNames()

const 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 PXDAnalyticGainCalibrationAlgorithm, PXDValidationAlgorithm, SVD3SampleCoGTimeCalibrationAlgorithm, SVD3SampleELSTimeCalibrationAlgorithm, SVDClusterAbsoluteTimeShifterAlgorithm, SVDCoGTimeCalibrationAlgorithm, TestBoundarySettingAlgorithm, and TestCalibrationAlgorithm.

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 503 of file CalibrationAlgorithm.cc.

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

◆ 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{
299 saveCalibration(data, name, m_data.getRequestedIov());
300}

◆ 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:150

◆ 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{
294 saveCalibration(data, name, m_data.getRequestedIov());
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 ( const 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 m_granularityOfData
Granularity of input data. This only changes when the input files change so it isn't specific to an e...
void setInputFileNames(PyObject *inputFileNames)
Set the input file names used for this algorithm from a Python list.
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.

◆ setInputFileNames() [2/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}
Scalar convertPythonObject(const boost::python::object &pyObject, Scalar)
Convert from Python to given type.

◆ setOutputJsonValue()

template<class T>
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,
const int run,
const int experiment )
staticprotectedinherited

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

Definition at line 405 of file CalibrationAlgorithm.cc.

406{
407 // Construct an EventMetaData object but NOT in the Datastore
408 EventMetaData emd(event, run, experiment);
409 // Explicitly update while avoiding registering a Datastore object
411 // Also update the intra-run objects to the event at the same time (maybe unnecessary...)
413}

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.

385{""};

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

388{""};

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