Belle II Software development
eclLeakageAlgorithm Class Reference

Calculate ECL energy leakage corrections. More...

#include <eclLeakageAlgorithm.h>

Inheritance diagram for eclLeakageAlgorithm:
CalibrationAlgorithm

Public Types

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

Public Member Functions

 eclLeakageAlgorithm ()
 Constructor.
 
virtual ~eclLeakageAlgorithm ()
 Destructor.
 
void setLowEnergyThreshold (double lowEnergyThreshold)
 Setter for m_lowEnergyThreshold.
 
double getLowEnergyThreshold ()
 Getter for m_lowEnergyThreshold.
 
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.
 
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

double m_lowEnergyThreshold = 0.0
 Parameters to control fit procedure.
 
int t_cellID = 0
 For TTree.
 
int t_thetaID = 0
 thetaID of photon
 
int t_region = 0
 region of photon 0=forward 1=barrel 2=backward
 
int t_thetaBin = -1
 binned location in theta relative to crystal edge
 
int t_phiBin = -1
 binned location in phi relative to crystal edge
 
int t_phiMech = -1
 mechanical structure next to lower phi (0), upper phi (1), or neither (2)
 
int t_energyBin = -1
 generated energy point
 
int t_nCrys = -1
 number of crystals used to calculate energy
 
float t_energyFrac = 0.
 measured energy after nOptimal bias and peak corrections, divided by generated
 
float t_origEnergyFrac = 0.
 corrected energy at time of generation, divided by generated
 
float t_locationError = 999.
 reconstructed minus generated position (cm)
 
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

Calculate ECL energy leakage corrections.

Definition at line 22 of file eclLeakageAlgorithm.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

◆ eclLeakageAlgorithm()

Constructor.


Definition at line 170 of file eclLeakageAlgorithm.cc.

170 : CalibrationAlgorithm("eclLeakageCollector")
171{
173 "Generate payload ECLLeakageCorrection from single photon MC recorded by eclLeakageCollectorModule"
174 );
175}
Base class for calibration algorithms.
void setDescription(const std::string &description)
Set algorithm description (in constructor)

◆ ~eclLeakageAlgorithm()

virtual ~eclLeakageAlgorithm ( )
inlinevirtual

Destructor.

Definition at line 29 of file eclLeakageAlgorithm.h.

29{}

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.


< number of thetaID

< 0 = forward, 1 = barrel, 2 = backward

< first thetaID of the barrel

< last thetaID of the barrel

< first thetaID to find leakage corrections

< last thetaID to find leakage correction

< locations across crystal

< number of energy points

< number of thetaID/energy combinations

< minimum allowed eta in Novosibirsk fits

< maximum allowed eta

< categories of fit status

< log(E) values for each region

Implements CalibrationAlgorithm.

Definition at line 180 of file eclLeakageAlgorithm.cc.

181{
182
183 //====================================================================================
184 //====================================================================================
185 //..Step 1. Set up prior to first loop through data
186 B2INFO("starting eclLeakageAlgorithm");
187
188 //-----------------------------------------------------------------------------------
189 //..Read in histograms created by the collector and fix normalization
190 auto inputParameters = getObjectPtr<TH1F>("inputParameters");
191 int lastBin = inputParameters->GetNbinsX();
192 double scale = inputParameters->GetBinContent(lastBin); // number of times inputParameters was filled
193 for (int ib = 1; ib < lastBin; ib++) {
194 double param = inputParameters->GetBinContent(ib);
195 inputParameters->SetBinContent(ib, param / scale);
196 inputParameters->SetBinError(ib, 0.);
197 }
198
199 //..And write to disk.
200 TFile* histfile = new TFile("eclLeakageAlgorithm.root", "recreate");
201 inputParameters->Write();
202
203 //..Parameters
204 const int nThetaID = 69;
205 const int nLeakReg = 3;
206 const int firstBarrelID = 13;
207 const int lastBarrelID = 58;
208 const int firstUsefulThID = 3;
209 const int lastUsefulThID = 66;
211 const int nPositions = (int)(inputParameters->GetBinContent(1) + 0.001);
212 const int nEnergies = (int)(inputParameters->GetBinContent(2) + 0.001);
213 const int nbinX = nEnergies * nThetaID;
215 const double etaMin = -5.;
216 const double etaMax = 5.;
219 //..Energies
220 auto generatedE = create2Dvector<float>(nLeakReg, nEnergies); // 3 regions forward barrel backward
221
222 int bin = 2; // bin 1 = nPositions, bin 2 = nEnergies, bin 3 = first energy
223 for (int ireg = 0; ireg < nLeakReg; ireg++) {
224 B2INFO("Generated energies for ireg = " << ireg);
225 for (int ie = 0; ie < nEnergies; ie++) {
226 bin++;
227 generatedE[ireg][ie] = inputParameters->GetBinContent(bin);
228 B2INFO(" " << ie << " " << generatedE[ireg][ie] << " GeV");
229 }
230 }
231 B2INFO("Low energy threshold = " << m_lowEnergyThreshold);
232
233 //..Energy per thetaID (in MeV, for use in titles etc)
234 auto iEnergiesMeV = create2Dvector<int>(nEnergies, nThetaID);
235 for (int thID = 0; thID < nThetaID; thID++) {
236 int ireg = 0;
237 if (thID >= firstBarrelID and thID <= lastBarrelID) {
238 ireg = 1;
239 } else if (thID > lastBarrelID) {
240 ireg = 2;
241 }
242 for (int ie = 0; ie < nEnergies; ie++) {
243 iEnergiesMeV[ie][thID] = (int)(1000.*generatedE[ireg][ie] + 0.5);
244 }
245 }
246
247 //-----------------------------------------------------------------------------------
248 //..Bins for eFrac histograms (eFrac = uncorrected reconstructed E / E true)
249 const double eFracLo = 0.4; // low edge of eFrac histograms
250 const double eFracHi = 1.5; // high edge of eFrac histograms
251 auto nEfracBins = create2Dvector<int>(nEnergies, nThetaID);
252 for (int thID = 0; thID < nThetaID; thID++) {
253 B2DEBUG(25, "eFrac nBins for thetaID " << thID);
254 for (int ie = 0; ie < nEnergies; ie++) {
255
256 //..ballpark resolution
257 double res_squared = 0.0001 + 0.064 / iEnergiesMeV[ie][thID];
258
259 //..Convert this to an even integer to get number of bins
260 double binNumOver2 = 3. / sqrt(res_squared);
261 int tempNBin = (int)(binNumOver2 + 0.5);
262 nEfracBins[ie][thID] = 2 * tempNBin;
263 B2DEBUG(25, " ie = " << ie << " E = " << iEnergiesMeV[ie][thID] << " nBins = " << nEfracBins[ie][thID]);
264 }
265 }
266
267 //-----------------------------------------------------------------------------------
268 //..Set up Novosibirsk fit function
269 TString statusString[7] = {"good", "refit", "lowStat", "lowProb", "peakAtLimit", "sigmaAtLimit", "etaAtLimit"};
270 const double fracEnt[2] = {0.683, 0.5}; // fit range includes 68% or 50% of entries
271 const double minEntries = 100.; // don't use fits with fewer entries
272 const double minMaxBin = 50.; // rebin if max bin is below this value
273 const double highMaxBin = 300.; // can float eta if max bin is above this value
274
275 TF1* func = new TF1("eclLeakageNovo", eclLeakageNovo, eFracLo, eFracHi, 4);
276 func->SetParNames("normalization", "peak", "effSigma", "eta");
277 func->SetLineColor(kRed);
278
279 //-----------------------------------------------------------------------------------
280 //..Specify the TTree
281 auto tree = getObjectPtr<TTree>("tree");
282 tree->SetBranchAddress("cellID", &t_cellID);
283 tree->SetBranchAddress("thetaID", &t_thetaID);
284 tree->SetBranchAddress("region", &t_region);
285 tree->SetBranchAddress("thetaBin", &t_thetaBin);
286 tree->SetBranchAddress("phiBin", &t_phiBin);
287 tree->SetBranchAddress("phiMech", &t_phiMech);
288 tree->SetBranchAddress("energyBin", &t_energyBin);
289 tree->SetBranchAddress("nCrys", &t_nCrys);
290 tree->SetBranchAddress("energyFrac", &t_energyFrac);
291 tree->SetBranchAddress("origEnergyFrac", &t_origEnergyFrac);
292 tree->SetBranchAddress("locationError", &t_locationError);
293
294 const int treeEntries = tree->GetEntries();
295 B2INFO("eclLeakageAlgorithm entries = " << treeEntries);
296
297 //-----------------------------------------------------------------------------------
298 //..Get current payload to help validate the new payload
299
300 //..Set event, run, exp number
301 const auto exprun = getRunList();
302 const int iEvt = 1;
303 const int iRun = exprun[0].second;
304 const int iExp = exprun[0].first;
307 evtPtr.registerInDataStore();
309 evtPtr.construct(iEvt, iRun, iExp);
310 DBStore& dbstore = DBStore::Instance();
311 dbstore.update();
312
313 //..Existing payload
314 DBObjPtr<ECLLeakageCorrections> existingCorrections("ECLLeakageCorrections");
315 TH2F existingThetaCorrection = existingCorrections->getThetaCorrections();
316 existingThetaCorrection.SetName("existingThetaCorrection");
317 TH2F existingPhiCorrection = existingCorrections->getPhiCorrections();
318 existingPhiCorrection.SetName("existingPhiCorrection");
319
320 //..Write out the correction histograms
321 histfile->cd();
322 existingThetaCorrection.Write();
323 existingPhiCorrection.Write();
324
325
326 //====================================================================================
327 //====================================================================================
328 //..Step 2. First loop, fill histograms of e/eTrue for each energy and thetaID
329
330 //-----------------------------------------------------------------------------------
331 //..One histogram of e/eTrue per energy per thetaID.
332 // Used to fix eta in subsequent fits, and to get overall correction for that thetaID
333 // (which should be very close to 1).
334 TString name;
335 TString title;
336 auto hELabUncorr = create2Dvector<TH1F*>(nEnergies, nThetaID);
337 for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
338 TString sthID = std::to_string(thID);
339 for (int ie = 0; ie < nEnergies; ie++) {
340 name = "hELabUncorr_" + std::to_string(ie) + "_" + sthID;
341 title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + ";E/Etrue";
342
343 //..High statistics for these plots; use more bins
344 hELabUncorr[ie][thID] = new TH1F(name, title, 2 * nEfracBins[ie][thID], eFracLo, eFracHi);
345 }
346 }
347
348 //-----------------------------------------------------------------------------------
349 //..Loop over events and store eFrac
350 for (int i = 0; i < treeEntries; i++) {
351 tree->GetEntry(i);
352
353 //..Only events with good reconstruction
354 bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
355 if (not goodReco) {continue;}
356
357 //..Fill histogram for full thetaID
358 hELabUncorr[t_energyBin][t_thetaID]->Fill(t_energyFrac);
359 }
360
361 //-----------------------------------------------------------------------------------
362 //..Fit each thetaID/energy histogram to get peak (overall correction) and eta (fixed
363 // in subsequent fits to individual locations).
364 auto peakUncorr = create2Dvector<float>(nEnergies, nThetaID); // store peak from each fit
365 auto etaUncorr = create2Dvector<float>(nEnergies, nThetaID); // store eta from each fit
366
367 std::vector<TString> failedELabUncorr; // names of hists with failed fits
368 std::vector<int> statusELabUncorr; // status of failed fits
369 int payloadStatus = 0; // Overall status of payload determination
370
371 //..Record fit status
372 TH1F* statusOfhELabUncorr = new TH1F("statusOfhELabUncorr",
373 "status of hELabUncorr fits for each thetaID/energy. 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb, 4 = peakAtLimit, 5 = sigmaAtLimit",
374 6, 0,
375 6);
376
377
378 //..Fit each histogram
379 for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
380 TString sthID = std::to_string(thID);
381 for (int ie = 0; ie < nEnergies; ie++) {
382 TH1F* hEnergy = (TH1F*)hELabUncorr[ie][thID];
383 double peak = -1.;
384 double eta = 0.;
385 int fitStatus = 2;
386 double entries = hEnergy->Integral();
387 int nIter = 0; // keep track of attempts to fit this histogram
388 double genE = iEnergiesMeV[ie][thID] / 1000.;
389 bool fitHist = entries > minEntries;
390
391 //..Possibly iterate fit starting from this point
392 while (fitHist) {
393
394 //..Fit parameters
395 double norm = hEnergy->GetMaximum();
396 double target = fracEnt[nIter] * entries; // fit range contains 68% or 50%
397 std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
398 startingParameters = eclLeakageFitParameters(hEnergy, target);
399 func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
400 func->SetParLimits(1, startingParameters[3], startingParameters[4]);
401 func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
402 func->SetParLimits(3, etaMin, etaMax);
403
404 //..Fit
405 name = hEnergy->GetName();
406 B2DEBUG(40, "Fitting " << name.Data());
407 hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
408 peak = func->GetParameter(1);
409 double effSigma = func->GetParameter(2);
410 eta = func->GetParameter(3);
411 double prob = func->GetProb();
412
413 //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
414 // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
415 double dEta = min((etaMax - eta), (eta - etaMin));
416 fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
417
418 //..If the fit probability is low, refit using a smaller range (fracEnt)
419 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
420 nIter++;
421 } else {
422 fitHist = false;
423 }
424 }
425
426 //-----------------------------------------------------------------------------------
427 //..Record failures and fit results.
428 // Mark payload as failed if energy is above low-energy threshold.
429 statusOfhELabUncorr->Fill(fitStatus + 0.000001);
430 if (fitStatus == 2 or fitStatus >= 4) {
431 statusELabUncorr.push_back(fitStatus);
432 failedELabUncorr.push_back(hEnergy->GetName());
433 if (genE > m_lowEnergyThreshold) {
434 payloadStatus = 1; // algorithm failed
435 } else {
436 peak = -1.; // fix up later
437 }
438 }
439 peakUncorr[ie][thID] = peak;
440 etaUncorr[ie][thID] = eta;
441
442 //..Write to disk
443 if (entries > minEntries) {
444 histfile->cd();
445 hELabUncorr[ie][thID]->Write();
446 }
447 }
448 }
449
450 //..Write out summary of fit status
451 histfile->cd();
452 statusOfhELabUncorr->Write();
453
454 //-----------------------------------------------------------------------------------
455 //..Quit now if one of the high-energy fits failed, since we will not get a payload.
456 int nbadFit = statusELabUncorr.size();
457 if (nbadFit > 0) {B2ERROR("hELabUncorr fit failures (one histogram per energy/thetaID):");}
458 for (int ibad = 0; ibad < nbadFit; ibad++) {
459 int badStat = statusELabUncorr[ibad];
460 B2ERROR(" histogram " << failedELabUncorr[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
461 }
462 if (payloadStatus != 0) {
463 B2ERROR("ecLeakageAlgorithm: fit to hELabUncorr failed. ");
464 histfile->Close();
465 return c_Failure;
466 }
467
468 //-----------------------------------------------------------------------------------
469 //..Fix up any failed (low-energy) fits by using a neighbouring thetaID
470 for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
471 for (int ie = 0; ie < nEnergies; ie++) {
472 if (peakUncorr[ie][thID] < 0.) {
473 if (thID > 40) {
474 for (int nextID = thID - 1; thID >= firstUsefulThID; thID--) {
475 if (peakUncorr[ie][nextID] > 0.) {
476 peakUncorr[ie][thID] = peakUncorr[ie][nextID];
477 break;
478 }
479 }
480 } else {
481 for (int nextID = thID + 1; thID <= lastUsefulThID; thID++) {
482 if (peakUncorr[ie][nextID] > 0.) {
483 peakUncorr[ie][thID] = peakUncorr[ie][nextID];
484 break;
485 }
486 }
487 }
488 }
489
490 //..If we were unable to get a successful fit from a neighbour, give up
491 if (peakUncorr[ie][thID] < 0.) {
492 B2ERROR("ecLeakageAlgorithm: unable to correct hELabUncorr for thetaID " << thID << " ie " << ie);
493 histfile->Close();
494 return c_Failure;
495 }
496 }
497 }
498
499 //====================================================================================
500 //====================================================================================
501 //..Step 3. Second loop, fill histograms of e/eTrue as a function of position.
502 // 29 locations in theta, 29 in phi.
503 // Crystals next to mechanical structure in phi are treated separately from
504 // crystals without mechanical structure.
505
506 //-----------------------------------------------------------------------------------
507 //..Histograms to store the energy
508 const int nDir = 3;
509 const TString dirName[nDir] = {"theta", "phiMech", "phiNoMech"};
510
511 auto eFracPosition = create4Dvector<TH1F*>(nEnergies, nThetaID, nDir, nPositions); // the histograms
512
513
514 std::vector<TString> failedeFracPosition; // names of hists with failed fits
515 std::vector<int> statuseFracPosition; // status of failed fits
516
517 for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
518 TString sthID = std::to_string(thID);
519 for (int ie = 0; ie < nEnergies; ie++) {
520 TString sie = std::to_string(ie);
521 for (int idir = 0; idir < nDir; idir++) {
522 TString sidir = std::to_string(idir);
523 for (int ipos = 0; ipos < nPositions; ipos++) {
524 TString sipos = std::to_string(ipos);
525 name = "eFracPosition_" + sie + "_" + sthID + "_" + sidir + "_" + sipos;
526 title = "eFrac " + to_string(iEnergiesMeV[ie][thID]) + " MeV thetaID " + sthID + " " + dirName[idir] + " pos " + sipos +
527 "; E/Etrue";
528 eFracPosition[ie][thID][idir][ipos] = new TH1F(name, title, nEfracBins[ie][thID], eFracLo, eFracHi);
529 }
530 }
531 }
532 }
533
534 //..And some summary histograms
535 TH1F* statusOfeFracPosition = new TH1F("statusOfeFracPosition",
536 "eFrac fit status: -2 noData, -1 lowE, 0 good, 1 redo fit, 2 lowStat, 3 lowProb, 4 peakAtLimit, 5 sigmaAtLimit", 8, -2,
537 6);
538 TH1F* probOfeFracPosition = new TH1F("probOfeFracPosition", "fit probability of eFrac fits for each position;probability", 100, 0,
539 1);
540 TH1F* maxOfeFracPosition = new TH1F("maxOfeFracPosition", "max entries of eFrac histograms;maximum bin content", 100, 0, 1000);
541
542 //-----------------------------------------------------------------------------------
543 //..Loop over events and store eFrac
544 for (int i = 0; i < treeEntries; i++) {
545 tree->GetEntry(i);
546
547 //..Only events with good reconstruction
548 bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
549 if (not goodReco) {continue;}
550
551 //..Theta location (idir = 0)
552 int idir = 0;
553 eFracPosition[t_energyBin][t_thetaID][idir][t_thetaBin]->Fill(t_energyFrac);
554
555 //..phi location. Note that t_phiBin starts at mechanical edge, if there is one.
556 idir = t_phiMech + 1;
557 eFracPosition[t_energyBin][t_thetaID][idir][t_phiBin]->Fill(t_energyFrac);
558 }
559
560 //-----------------------------------------------------------------------------------
561 //..Now fit many many histograms to get the position-dependent corrections
562 auto positionCorrection = create4Dvector<float>(nEnergies, nThetaID, nDir, nPositions);
563 auto positionCorrectionUnc = create4Dvector<float>(nEnergies, nThetaID, nDir, nPositions);
564 int nHistToFit = 0;
565
566
567 //..Temp histogram of position corrections
568 TH1F* thetaCorSummary = new TH1F("thetaCorSummary", "Theta dependent corrections;theta dependent correction", 100, 0.4, 1.4);
569 TH1F* phiCorSummary = new TH1F("phiCorSummary", "Phi dependent corrections;phi dependent correction", 100, 0.4, 1.4);
570
571 for (int thID = firstUsefulThID; thID <= lastUsefulThID; thID++) {
572 for (int ie = 0; ie < nEnergies; ie++) {
573 double genE = iEnergiesMeV[ie][thID] / 1000.;
574 for (int idir = 0; idir < nDir; idir++) {
575 for (int ipos = 0; ipos < nPositions; ipos++) {
576 TH1F* hEnergy = (TH1F*)eFracPosition[ie][thID][idir][ipos];
577 if (hEnergy->Integral() > minEntries) {nHistToFit++;}
578
579 //..Default peak from fit to full thetaID/energy = peakUncorr;
580 // correction = peak / sqrt(peakUncorr) = sqrt(peakUncorr)
581 double correction = sqrt(peakUncorr[ie][thID]);
582 double corrUnc = 0.05; // arbitrary uncertainty
583 double prob = -1.;
584 int fitStatus = 2; // low stats
585 if (genE < m_lowEnergyThreshold) {fitStatus = -1;} // low energy, don't fit, use default peak value
586 if (hEnergy->GetEntries() < 0.5) {fitStatus = -2;} // unused, eg barrel pos=2
587 int nIter = 0; // keep track of attempts to fit this histogram
588 double entries = hEnergy->Integral();
589 bool fitHist = entries > minEntries and genE >= m_lowEnergyThreshold;
590
591 //..Possibly iterate fit starting from this point
592 while (fitHist) {
593 fitHist = false;
594
595 //..Fit parameters
596 double target = fracEnt[nIter] * entries; // fit range contains 68%
597 std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
598 bool findParm = true;
599 while (findParm) {
600 startingParameters = eclLeakageFitParameters(hEnergy, target);
601
602 //..Rebin if stats are low and we have enough DOF
603 if (hEnergy->GetMaximum()<minMaxBin and startingParameters[5]>11.9) {
604 hEnergy->Rebin(2);
605 } else {
606 findParm = false;
607 }
608 }
609 double norm = hEnergy->GetMaximum(); // max bin after rebinning
610 double fitLow = startingParameters[3];
611 double fitHigh = startingParameters[4];
612
613 //..Eta from the fit to the full crystal
614 double etaFix = etaUncorr[ie][thID];
615 func->SetParameters(norm, startingParameters[0], startingParameters[1], etaFix);
616 func->SetParLimits(1, fitLow, fitHigh);
617 func->SetParLimits(2, 0., fitHigh - fitLow);
618
619 //..Fix eta, unless really good statistics
620 if (hEnergy->GetMaximum() < highMaxBin) {
621 func->SetParLimits(3, etaFix, etaFix);
622 } else {
623 func->SetParLimits(3, etaMin, etaMax);
624 }
625
626 //..Fit
627 name = hEnergy->GetName();
628 B2DEBUG(40, "Fitting " << name.Data());
629 hEnergy->Fit(func, "LIEQ", "", fitLow, fitHigh);
630 double peak = func->GetParameter(1);
631 double effSigma = func->GetParameter(2);
632 double eta = func->GetParameter(3);
633 prob = func->GetProb();
634
635 //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
636 // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit.
637 double dEta = min((etaMax - eta), (eta - etaMin));
638 fitStatus = eclLeakageFitQuality(fitLow, fitHigh, peak, effSigma, dEta, prob);
639
640 //..If the fit probability is low, refit using a smaller range (fracEnt)
641 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
642 nIter++;
643 fitHist = true;
644
645 // Store correction except for peak or sigma at limit.
646 } else if (fitStatus <= 3) {
647
648 //..Divide by sqrt(peak for full crystal) to get correct average
649 // when multiplying the theta and phi corrections
650 correction = peak / sqrt(peakUncorr[ie][thID]);
651 corrUnc = func->GetParError(1) / sqrt(peakUncorr[ie][thID]);
652 }
653 }
654
655 //..Store the correction for this position
656 positionCorrection[ie][thID][idir][ipos] = correction;
657 positionCorrectionUnc[ie][thID][idir][ipos] = corrUnc;
658 statusOfeFracPosition->Fill(fitStatus + 0.00001);
659 probOfeFracPosition->Fill(prob);
660 maxOfeFracPosition->Fill(hEnergy->GetMaximum());
661
662 //..Summary histogram for successful fits. Note that fitStatus<0 also has prob<0.
663 if (prob > 0 and fitStatus <= 3) {
664 if (idir == 0) {
665 thetaCorSummary->Fill(correction);
666 } else {
667 phiCorSummary->Fill(correction);
668 }
669 }
670
671 //..Record the failed fit details. We may still use the correction.
672 if (fitStatus >= 2) {
673 statuseFracPosition.push_back(fitStatus);
674 failedeFracPosition.push_back(hEnergy->GetName());
675 histfile->cd();
676 hEnergy->Write();
677 }
678 }
679 }
680 }
681 }
682
683 //-----------------------------------------------------------------------------------
684 //..Summarize position fit results
685 B2INFO("eclLeakageAlgorithm: " << nHistToFit << " eFracPosition histograms to fit");
686 nbadFit = statuseFracPosition.size();
687 B2INFO("eFracPosition failed fits: " << nbadFit);
688 for (int ibad = 0; ibad < nbadFit; ibad++) {
689 int badStat = statuseFracPosition[ibad];
690 B2ERROR(" histogram " << failedeFracPosition[ibad].Data() << " status " << badStat << " " << statusString[badStat].Data());
691 }
692
693 //..Write to disk
694 histfile->cd();
695 statusOfeFracPosition->Write();
696 probOfeFracPosition->Write();
697 maxOfeFracPosition->Write();
698 thetaCorSummary->Write();
699 phiCorSummary->Write();
700
701
702 //====================================================================================
703 //====================================================================================
704 //..Step 4. Store quantities needed for the payloads
705
706 //-----------------------------------------------------------------------------------
707 //..First, we need to fix up the thetaID that have been so far missed.
708 // Just copy the values from the first or last thetaID for which corrections were found.
709
710 //..ThetaID before the first useful one
711 for (int thID = firstUsefulThID - 1; thID >= 0; thID--) {
712 for (int ie = 0; ie < nEnergies; ie++) {
713
714 for (int idir = 0; idir < 3; idir++) {
715 for (int ipos = 0; ipos < nPositions; ipos++) {
716 positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID + 1][idir][ipos];
717 positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID + 1][idir][ipos];
718 }
719 }
720 }
721 }
722
723 //..ThetaID beyond last useful one
724 for (int thID = lastUsefulThID + 1; thID < nThetaID; thID++) {
725 for (int ie = 0; ie < nEnergies; ie++) {
726
727 for (int idir = 0; idir < 3; idir++) {
728 for (int ipos = 0; ipos < nPositions; ipos++) {
729 positionCorrection[ie][thID][idir][ipos] = positionCorrection[ie][thID - 1][idir][ipos];
730 positionCorrectionUnc[ie][thID][idir][ipos] = positionCorrectionUnc[ie][thID - 1][idir][ipos];
731 }
732 }
733 }
734 }
735
736 //-----------------------------------------------------------------------------------
737 //..std::vectors of generated energies
738 std::vector<float> forwardVector;
739 std::vector<float> barrelVector;
740 std::vector<float> backwardVector;
741 for (int ie = 0; ie < nEnergies; ie++) {
742 forwardVector.push_back(log(generatedE[0][ie]));
743 barrelVector.push_back(log(generatedE[1][ie]));
744 backwardVector.push_back(log(generatedE[2][ie]));
745 }
746
747 //..Store in array format for validation studies
748 float leakLogE[3][12] = {};
749 for (int ie = 0; ie < nEnergies; ie++) {
750 leakLogE[0][ie] = forwardVector[ie];
751 leakLogE[1][ie] = barrelVector[ie];
752 leakLogE[2][ie] = backwardVector[ie];
753 }
754
755 //-----------------------------------------------------------------------------------
756 //..Position dependent corrections
757
758 //..2D histogram of theta-dependent corrections
759 TH2F thetaCorrection("thetaCorrection", "Theta correction vs bin;bin = thetaID + 69*energyBin;local theta bin", nbinX, 0, nbinX,
760 nPositions, 0, nPositions);
761 int ix = 0;
762 for (int ie = 0; ie < nEnergies; ie++) {
763 for (int thID = 0; thID < nThetaID; thID++) {
764 ix++;
765 for (int ipos = 0; ipos < nPositions; ipos++) {
766 int iy = ipos + 1;
767 float correction = positionCorrection[ie][thID][0][ipos];
768 float corrUnc = positionCorrectionUnc[ie][thID][0][ipos];
769 thetaCorrection.SetBinContent(ix, iy, correction);
770 thetaCorrection.SetBinError(ix, iy, corrUnc);
771 }
772 }
773 }
774
775 //..2D histogram of phi-dependent corrections. Twice as many x bins;
776 // one set for crystals next to mechanical structure, one for otherwise.
777 TH2F phiCorrection("phiCorrection", "Phi correction vs bin;bin = thetaID + 69*energyBin;local phi bin", 2 * nbinX, 0, 2 * nbinX,
778 nPositions, 0, nPositions);
779 ix = 0;
780 for (int ie = 0; ie < nEnergies; ie++) {
781 for (int thID = 0; thID < nThetaID; thID++) {
782 ix++;
783 for (int ipos = 0; ipos < nPositions; ipos++) {
784 int iy = ipos + 1;
785
786 //..First set of corrections
787 float correction = positionCorrection[ie][thID][1][ipos];
788 float corrUnc = positionCorrectionUnc[ie][thID][1][ipos];
789 phiCorrection.SetBinContent(ix, iy, correction);
790 phiCorrection.SetBinError(ix, iy, corrUnc);
791
792 //..Second set
793 correction = positionCorrection[ie][thID][2][ipos];
794 corrUnc = positionCorrectionUnc[ie][thID][2][ipos];
795 phiCorrection.SetBinContent(ix + nbinX, iy, correction);
796 phiCorrection.SetBinError(ix + nbinX, iy, corrUnc);
797 }
798 }
799 }
800
801 //-----------------------------------------------------------------------------------
802 //..nCrys dependent corrections; not used since nOptimal was implemented.
803 const int maxN = 21; // maximum number of crystals in a cluster
804 TH2F nCrystalCorrection("nCrystalCorrection", "nCrys correction vs bin;bin = thetaID + 69*energyBin;nCrys", nbinX, 0, nbinX,
805 maxN + 1, 0, maxN + 1);
806
807 ix = 0;
808 for (int ie = 0; ie < nEnergies; ie++) {
809 for (int thID = 0; thID < nThetaID; thID++) {
810 ix++;
811 for (int in = 0; in <= maxN; in++) {
812 int iy = in + 1;
813 float correction = 1.;
814 float corrUnc = 0.;
815 nCrystalCorrection.SetBinContent(ix, iy, correction);
816 nCrystalCorrection.SetBinError(ix, iy, corrUnc);
817 }
818 }
819 }
820
821
822 //====================================================================================
823 //====================================================================================
824 //..Step 5. Diagnostic histograms
825
826 //..Start by storing the payload histograms
827 histfile->cd();
828 thetaCorrection.Write();
829 phiCorrection.Write();
830 nCrystalCorrection.Write();
831
832 //-----------------------------------------------------------------------------------
833 //..One histogram of new and original reconstructed energy after leakage correction
834 // per generated energy per region. Also uncorrected.
835 // The "Corrected no nCrys" and "Corrected measured" are identical, but keep both
836 // for easier comparison with earlier results.
837 const int nResType = 5;
838 const TString resName[nResType] = {"Uncorrected", "Original", "Corrected no nCrys", "Corrected measured", "Corrected true"};
839 const TString regName[nLeakReg] = {"forward", "barrel", "backward"};
840 auto energyResolution = create3Dvector<TH1F*>(nLeakReg, nEnergies, nResType);
841
842 //..Base number of bins on a typical thetaID for each region
843 int thIDReg[nLeakReg];
844 thIDReg[0] = 9;
845 thIDReg[1] = 35;
846 thIDReg[2] = 61;
847
848 for (int ireg = 0; ireg < nLeakReg; ireg++) {
849 for (int ie = 0; ie < nEnergies; ie++) {
850
851 //..Titles and bins for this region and energy
852 TString stireg = std::to_string(ireg);
853 TString stie = std::to_string(ie);
854 TString stieName = std::to_string(iEnergiesMeV[ie][thIDReg[ireg]]);
855 int nbinReg = 20 * nEfracBins[ie][thIDReg[ireg]];
856
857 for (int ires = 0; ires < nResType; ires++) {
858 name = "energyResolution_" + stireg + "_" + stie + "_" + std::to_string(ires);
859 title = resName[ires] + " energy, " + stieName + " MeV, " + regName[ireg] + ";Reconstructed E/Etrue";
860 energyResolution[ireg][ie][ires] = new TH1F(name, title, nbinReg, eFracLo, eFracHi);
861 }
862 }
863 }
864
865 //-----------------------------------------------------------------------------------
866 //..Loop over events and store leakage-corrected energies
867 for (int i = 0; i < treeEntries; i++) {
868 tree->GetEntry(i);
869
870 //..Only events with good reconstruction
871 bool goodReco = t_thetaID >= firstUsefulThID and t_thetaID <= lastUsefulThID and t_energyBin >= 0;
872 if (not goodReco) {continue;}
873
874 //-----------------------------------------------------------------------------------
875 //..Corrections using true energy
876
877 //..Position-dependent leakage corrections using true energy
878 int idir = 0;
879 float thetaPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_thetaBin];
880
881 idir = t_phiMech + 1;
882 float phiPosCor = positionCorrection[t_energyBin][t_thetaID][idir][t_phiBin];
883 float corrTrue = thetaPosCor * phiPosCor;
884
885 //-----------------------------------------------------------------------------------
886 //..Find correction using measured energy. The correction is a function of corrected
887 // energy, so will need to iterate
888 float corrMeasured = 0.96; // typical correction as starting point
889 for (int iter = 0; iter < 2; iter++) {
890
891
892 //..Energy points that bracket this value
893 float energyRaw = t_energyFrac * generatedE[t_region][t_energyBin];
894 float logEnergy = log(energyRaw / corrMeasured);
895 int ie0 = 0; // lower energy point
896 if (logEnergy < leakLogE[t_region][0]) {
897 ie0 = 0;
898 } else if (logEnergy > leakLogE[t_region][nEnergies - 1]) {
899 ie0 = nEnergies - 2;
900 } else {
901 while (logEnergy > leakLogE[t_region][ie0 + 1]) {ie0++;}
902 }
903
904 //..Corrections from lower and upper energy points
905 float cor0 = positionCorrection[ie0][t_thetaID][0][t_thetaBin] * positionCorrection[ie0][t_thetaID][t_phiMech + 1][t_phiBin];
906 float cor1 = positionCorrection[ie0 + 1][t_thetaID][0][t_thetaBin] * positionCorrection[ie0 + 1][t_thetaID][t_phiMech +
907 1][t_phiBin];
908
909 //..Interpolate in logE
910 corrMeasured = cor0 + (cor1 - cor0) * (logEnergy - leakLogE[t_region][ie0]) / (leakLogE[t_region][ie0 + 1] -
911 leakLogE[t_region][ie0]);
912 }
913
914 //..No longer have a separate case that excludes the nCrys correction
915 float corrNonCrys = corrMeasured;
916
917 //-----------------------------------------------------------------------------------
918 //..Fill the histograms
919 energyResolution[t_region][t_energyBin][0]->Fill(t_energyFrac); // uncorrected
920 energyResolution[t_region][t_energyBin][1]->Fill(t_origEnergyFrac); // original payload
921 energyResolution[t_region][t_energyBin][2]->Fill(t_energyFrac / corrNonCrys); // no nCrys, measured energy
922 energyResolution[t_region][t_energyBin][3]->Fill(t_energyFrac / corrMeasured); // corrected, measured energy
923 energyResolution[t_region][t_energyBin][4]->Fill(t_energyFrac / corrTrue); // corrected, true energy
924 }
925
926 //-----------------------------------------------------------------------------------
927 //..Fit each histogram to find peak, and extract resolution.
928
929 //..Store the peak and resolution values for each histogram
930 auto peakEnergy = create3Dvector<float>(nLeakReg, nEnergies, nResType);
931 auto energyRes = create3Dvector<float>(nLeakReg, nEnergies, nResType);
932
933 //..Loop over the histograms to be fit
934 for (int ireg = 0; ireg < nLeakReg; ireg++) {
935 for (int ie = 0; ie < nEnergies; ie++) {
936
937 //..Base the resolution on the number of entries in the corrected plot for all 4 cases.
938 double entries = energyResolution[ireg][ie][nResType - 1]->Integral();
939 for (int ires = 0; ires < nResType; ires++) {
940 TH1F* hEnergy = (TH1F*)energyResolution[ireg][ie][ires];
941 double peak = -1.;
942 double res68 = 99.; // resolution based on 68.3% of the entries
943 int nIter = 0; // keep track of attempts to fit this histogram
944 bool fitHist = entries > minEntries;
945
946 //..Possibly iterate fit starting from this point
947 while (fitHist) {
948
949 //..Fit parameters
950 double norm = hEnergy->GetMaximum();
951 double target = fracEnt[nIter] * entries; // fit range contains 68.3% or 50%
952 std::vector<double> startingParameters;// peak, sigma, eta, fitLow, fitHigh, nbins
953 startingParameters = eclLeakageFitParameters(hEnergy, target);
954 func->SetParameters(norm, startingParameters[0], startingParameters[1], startingParameters[2]);
955 func->SetParLimits(1, startingParameters[3], startingParameters[4]);
956 func->SetParLimits(2, 0., startingParameters[4] - startingParameters[3]);
957 func->SetParLimits(3, etaMin, etaMax);
958
959 //..First iteration the fitting range contains >68.3% of the histogram integral.
960 if (nIter == 0) {
961
962 //..Find the bin numbers used in the fit
963 int minLo = hEnergy->GetXaxis()->FindBin(startingParameters[3] + 0.0001);
964 int minHi = hEnergy->GetXaxis()->FindBin(startingParameters[4] - 0.0001);
965
966 //..Subtract a partial bin to get to exactly 68.3%
967 double dx = hEnergy->GetBinLowEdge(minLo + 1) - hEnergy->GetBinLowEdge(minLo);
968 double overage = hEnergy->Integral(minLo, minHi) - target;
969 double subLo = overage / hEnergy->GetBinContent(minLo);
970 double subHi = overage / hEnergy->GetBinContent(minHi);
971
972 //..Resolution is half the range that contains 68.3% of the events
973 res68 = 0.5 * dx * (1 + minHi - minLo - max(subLo, subHi));
974 }
975
976 //..Fit
977 name = hEnergy->GetName();
978 B2DEBUG(40, "Fitting " << name.Data());
979 hEnergy->Fit(func, "LIEQ", "", startingParameters[3], startingParameters[4]);
980 peak = func->GetParameter(1);
981 double effSigma = func->GetParameter(2);
982 double eta = func->GetParameter(3);
983 double prob = func->GetProb();
984
985 //..Check fit quality 0 = good, 1 = redo fit, 2 = lowStat, 3 = lowProb,
986 // 4 = peakAtLimit, 5 = sigmaAtLimit, 6 = etaAtLimit
987 double dEta = min((etaMax - eta), (eta - etaMin));
988 int fitStatus = eclLeakageFitQuality(startingParameters[3], startingParameters[4], peak, effSigma, dEta, prob);
989
990 //..If the fit probability is low, refit using a smaller range (fracEnt)
991 if ((fitStatus == 1 or fitStatus == 3) and nIter == 0) {
992 nIter++;
993 } else {
994 fitHist = false;
995 }
996
997 } // end of while
998
999 //..Record peak and resolution
1000 peakEnergy[ireg][ie][ires] = peak;
1001 energyRes[ireg][ie][ires] = res68;
1002
1003 //..And write to disk
1004 histfile->cd();
1005 hEnergy->Write();
1006 }
1007 }
1008 }
1009
1010 //-----------------------------------------------------------------------------------
1011 //..Summarize resolution
1012 int nresBins = nEnergies * nLeakReg * (nResType + 1); // +1 to add an empty bin after each set
1013 TH1F* peakSummary = new TH1F("peakSummary", "Peak E/Etrue for each method, region, energy;Energy energy point", nresBins, 0,
1014 nEnergies);
1015 TH1F* resolutionSummary = new TH1F("resolutionSummary", "Resolution/peak for each method, region, energy;Energy energy point",
1016 nresBins, 0, nEnergies);
1017
1018 B2INFO("Resolution divided by peak for each energy bin and region " << nResType << " ways");
1019 for (int ires = 0; ires < nResType; ires++) {B2INFO(" " << resName[ires]);}
1020 ix = 0;
1021 for (int ie = 0; ie < nEnergies; ie++) {
1022 B2INFO(" energy point " << ie);
1023 for (int ireg = 0; ireg < nLeakReg; ireg++) {
1024 B2INFO(" region " << ireg);
1025 for (int ires = 0; ires < nResType; ires++) {
1026
1027 //..Store in summary histograms
1028 ix++;
1029 peakSummary->SetBinContent(ix, peakEnergy[ireg][ie][ires]);
1030 peakSummary->SetBinError(ix, 0.);
1031 resolutionSummary->SetBinContent(ix, energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1032 resolutionSummary->SetBinError(ix, 0.);
1033
1034 //..Print out as well
1035 B2INFO(" " << ires << " " << energyRes[ireg][ie][ires] / peakEnergy[ireg][ie][ires]);
1036 }
1037 ix++;
1038 }
1039 }
1040
1041 //..Write the summary histograms to disk
1042 histfile->cd();
1043 peakSummary->Write();
1044 resolutionSummary->Write();
1045
1046
1047 //====================================================================================
1048 //====================================================================================
1049 //..Step 6. Finish up.
1050
1051 //-----------------------------------------------------------------------------------
1052 //..Close histogram file
1053 histfile->Close();
1054
1055 //------------------------------------------------------------------------
1056 //..Create and store payload
1057 ECLLeakageCorrections* leakagePayload = new ECLLeakageCorrections();
1058 leakagePayload->setlogEnergiesFwd(forwardVector);
1059 leakagePayload->setlogEnergiesBrl(barrelVector);
1060 leakagePayload->setlogEnergiesBwd(backwardVector);
1061 leakagePayload->setThetaCorrections(thetaCorrection);
1062 leakagePayload->setPhiCorrections(phiCorrection);
1063 leakagePayload->setnCrystalCorrections(nCrystalCorrection);
1064 saveCalibration(leakagePayload, "ECLLeakageCorrections");
1065
1066 B2INFO("eclLeakageAlgorithm: successfully stored payload ECLLeakageCorrections");
1067 return c_OK;
1068}
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 leakage corrections, including nCrys dependence
void setThetaCorrections(const TH2F &thetaCorrections)
Set the 2D histogram containing the theta corrections for each thetaID and energy.
void setlogEnergiesBrl(const std::vector< float > &logEnergiesBrl)
Set the vector of energies used to evaluate the leakage corrections in the barrel.
void setnCrystalCorrections(const TH2F &nCrystalCorrections)
Set the 2D histogram containing the nCrys corrections for each thetaID and energy.
void setlogEnergiesFwd(const std::vector< float > &logEnergiesFwd)
Set the vector of energies used to evaluate the leakage corrections in the forward endcap.
void setPhiCorrections(const TH2F &phiCorrections)
Set the 2D histogram containing the phi corrections for each thetaID and energy.
void setlogEnergiesBwd(const std::vector< float > &logEnergiesBwd)
Set the vector of energies used to evaluate the leakage corrections in the backward endcap.
float t_energyFrac
measured energy after nOptimal bias and peak corrections, divided by generated
float t_locationError
reconstructed minus generated position (cm)
int t_energyBin
generated energy point
int t_phiBin
binned location in phi relative to crystal edge
int t_phiMech
mechanical structure next to lower phi (0), upper phi (1), or neither (2)
int t_region
region of photon 0=forward 1=barrel 2=backward
int t_thetaBin
binned location in theta relative to crystal edge
int t_nCrys
number of crystals used to calculate energy
double m_lowEnergyThreshold
Parameters to control fit procedure.
float t_origEnergyFrac
corrected energy at time of generation, divided by generated
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 update()
Updates all objects that are outside their interval of validity.
Definition: DBStore.cc:79
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

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

◆ getLowEnergyThreshold()

double getLowEnergyThreshold ( )
inline

Getter for m_lowEnergyThreshold.

Definition at line 35 of file eclLeakageAlgorithm.h.

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

◆ setLowEnergyThreshold()

void setLowEnergyThreshold ( double  lowEnergyThreshold)
inline

Setter for m_lowEnergyThreshold.

Definition at line 32 of file eclLeakageAlgorithm.h.

32{m_lowEnergyThreshold = lowEnergyThreshold;}

◆ 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
void updateEvent()
Updates all intra-run dependent objects.
Definition: DBStore.cc:142

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_lowEnergyThreshold

double m_lowEnergyThreshold = 0.0
private

Parameters to control fit procedure.

only minimal fits below this value

Definition at line 45 of file eclLeakageAlgorithm.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.

◆ t_cellID

int t_cellID = 0
private

For TTree.

cellID of photon

Definition at line 48 of file eclLeakageAlgorithm.h.

◆ t_energyBin

int t_energyBin = -1
private

generated energy point

Definition at line 54 of file eclLeakageAlgorithm.h.

◆ t_energyFrac

float t_energyFrac = 0.
private

measured energy after nOptimal bias and peak corrections, divided by generated

Definition at line 56 of file eclLeakageAlgorithm.h.

◆ t_locationError

float t_locationError = 999.
private

reconstructed minus generated position (cm)

Definition at line 58 of file eclLeakageAlgorithm.h.

◆ t_nCrys

int t_nCrys = -1
private

number of crystals used to calculate energy

Definition at line 55 of file eclLeakageAlgorithm.h.

◆ t_origEnergyFrac

float t_origEnergyFrac = 0.
private

corrected energy at time of generation, divided by generated

Definition at line 57 of file eclLeakageAlgorithm.h.

◆ t_phiBin

int t_phiBin = -1
private

binned location in phi relative to crystal edge

Definition at line 52 of file eclLeakageAlgorithm.h.

◆ t_phiMech

int t_phiMech = -1
private

mechanical structure next to lower phi (0), upper phi (1), or neither (2)

Definition at line 53 of file eclLeakageAlgorithm.h.

◆ t_region

int t_region = 0
private

region of photon 0=forward 1=barrel 2=backward

Definition at line 50 of file eclLeakageAlgorithm.h.

◆ t_thetaBin

int t_thetaBin = -1
private

binned location in theta relative to crystal edge

Definition at line 51 of file eclLeakageAlgorithm.h.

◆ t_thetaID

int t_thetaID = 0
private

thetaID of photon

Definition at line 49 of file eclLeakageAlgorithm.h.


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