Belle II Software  release-05-02-19
eclMuMuEAlgorithm Class Reference

Calibrate ecl crystals using muon pair events. More...

#include <eclMuMuEAlgorithm.h>

Inheritance diagram for eclMuMuEAlgorithm:
Collaboration diagram for eclMuMuEAlgorithm:

Public Types

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

Public Member Functions

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

Public Attributes

int cellIDLo
 ..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events More...
 
int cellIDHi
 Last cellID to be fit.
 
int minEntries
 All crystals to be fit must have at least minEntries events in the fit range.
 
int maxIterations
 no more than maxIteration iterations
 
double tRatioMin
 Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
 
double tRatioMax
 Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
 
bool performFits
 if false, input histograms are copied to output, but no fits are done.
 
bool findExpValues
 if true, fits are used to find expected energy deposit for each crystal instead of the calibration constant
 
int storeConst
 controls which values are written to the database. More...
 

Protected Member Functions

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

Protected Attributes

std::vector< Calibration::ExpRun > m_boundaries
 

Private Member Functions

std::string getExpRunString (Calibration::ExpRun &expRun) const
 Gets the "exp.run" string repr. of (exp,run)
 
std::string getFullObjectPath (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

int fitOK = 16
 fit is OK
 
int iterations = 8
 fit reached max number of iterations, but is useable
 
int atLimit = 4
 a parameter is at the limit; fit not useable
 
int poorFit = 3
 low chi square; fit not useable
 
int noPeak = 2
 Novosibirsk component of fit is negligible; fit not useable.
 
int notFit = -1
 no fit performed
 
std::vector< std::string > m_inputFileNames
 List of input files to the Algorithm, will initially be user defined but then gets the wildcards expanded during execute()
 
std::map< Calibration::ExpRun, std::vector< std::string > > m_runsToInputFiles
 Map of Runs to input files. Gets filled when you call getRunRangeFromAllData, gets cleared when setting input files again.
 
std::string m_granularityOfData
 Granularity of input data. This only changes when the input files change so it isn't specific to an execution.
 
ExecutionData m_data
 Data specific to a SINGLE execution of the algorithm. Gets reset at the beginning of execution.
 
std::string m_description {""}
 Description of the algorithm.
 
std::string m_prefix {""}
 The name of the TDirectory the collector objects are contained within.
 
nlohmann::json m_jsonExecutionInput = nlohmann::json::object()
 Optional input JSON object used to make decisions about how to execute the algorithm code.
 
nlohmann::json m_jsonExecutionOutput = nlohmann::json::object()
 Optional output JSON object that can be set during the execution by the underlying algorithm code.
 

Static Private Attributes

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

Detailed Description

Calibrate ecl crystals using muon pair events.

Definition at line 39 of file eclMuMuEAlgorithm.h.

Member Enumeration Documentation

◆ EResult

enum EResult
inherited

The result of calibration.

Enumerator
c_OK 

Finished successfuly =0 in Python.

c_Iterate 

Needs iteration =1 in Python.

c_NotEnoughData 

Needs more data =2 in Python.

c_Failure 

Failed =3 in Python.

c_Undefined 

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

Definition at line 50 of file CalibrationAlgorithm.h.

Member Function Documentation

◆ calibrate()

CalibrationAlgorithm::EResult calibrate ( )
overrideprotectedvirtual

..Run algorithm on events


ranges of various fit parameters, and tolerance to determine that fit is at the limit

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


Clean up existing histograms if necessary


Histograms containing the data collected by eclMuMuECollectorModule


Record the number of entries per crystal in the normalized energy histogram and calculate the average expected energy per crystal and calibration constants from Collector


Write out the basic histograms in all cases


If we have not been asked to do fits, we can quit now


Check that every crystal has enough entries.

Insufficient data. Quit if we are required to have a successful fit for every crystal


Some prep for the many fits about to follow

histograms to store results for database

Diagnostic histograms

1D summary histograms


Fits are requested and there is sufficient data. Loop over specified crystals and performs fits to the amplitude distributions

Project 1D histogram of energy in this crystal

Fit function (xmin, xmax, nparameters) for this histogram

Estimate initial parameters from the histogram. For peak, use maximum bin in the allowed range

Fit range is histogram low edge plus a few bins to peak + 2.5*effective sigma

eta and constant are nominal values

Constant from lower edge of plot

Parameters to control iterations. dIter checks if we are stuck in a loop


Iterate from this point

Set the initial parameters

Fix eta if this is an endcap crystal, or if findExpValues==false

Fit

The upper fit range should correspond to 20-25% of the peak. Iterate if necessary.

Check if we are oscillating between two end points

Many iterations may mean we are stuck in a loop. Try a different end point.

Set the constant term to 0 if we are close to the limit

No more than specified number of iterations


Fit status

No peak; normalization of Novo component is too small

poor fit, or relatively poor fit with too many iterations

parameter at limit

fill diagnostic histograms

1D summary histograms

Store histogram with fit


Interpret results of fit as expected energy or calibration constant

if the fit is not successful, set peakE to -1, so that output calib = -1 * abs(input calib)

Find expected energies from MC, if requested

Otherwise, calibration constant


Write output to DB if requested and successful

Store expected energies

Store calibration constants


Write out diagnostic histograms

Histograms containing values written to DB


Clean up histograms in case Algorithm is called again


Set the return code appropriately

Implements CalibrationAlgorithm.

Definition at line 51 of file eclMuMuEAlgorithm.cc.

52 {
55  double limitTol = 0.0005; /*< tolerance for checking if a parameter is at the limit */
56  double minFitLimit = 1e-25; /*< cut off for labeling a fit as poor */
57  double minFitProbIter = 1e-8; /*< cut off for labeling a fit as poor if it also has many iterations */
58  double constRatio = 0.5; /*< Novosibirsk normalization must be greater than constRatio x constant term */
59  double peakMin(0.4), peakMax(2.2); /*< range for peak of normalized energy distribution */
60  double peakTol = limitTol * (peakMax - peakMin); /*< fit is at limit if it is within peakTol of min or max */
61  double effSigMin(0.02), effSigMax(0.2); /*< range for effective sigma of normalized energy distribution */
62  double effSigTol = limitTol * (effSigMax - effSigMin); /*< fit is at limit if it is within effSigTol of min or max */
63  double etaNom(-0.41); /*< Nominal tail parameter; fixed to this value for endcap and when findExValues=false */
64  double etaMin(-0.7), etaMax(-0.2); /*< Novosibirsk tail parameter range */
65  double etaTol = limitTol * (etaMax - etaMin); /*< fit is at limit if it is within etaTol of min or max */
66  double constMin(0.), constMax(10.); /*< constant term in normalized energy distribution */
67  double constTol = limitTol * constMax; /*< if constant is less than constTol, it will be fixed to 0 */
68 
70  gROOT->SetBatch();
71 
74  TH1F* dummy;
75  dummy = (TH1F*)gROOT->FindObject("IntegralVsCrysID");
76  if (dummy) {delete dummy;}
77  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
78  if (dummy) {delete dummy;}
79  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
80  if (dummy) {delete dummy;}
81  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
82  if (dummy) {delete dummy;}
83 
86  auto TrkPerCrysID = getObjectPtr<TH1F>("TrkPerCrysID");
87  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
88  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
89  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
90  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
91  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
92  auto RawDigitAmpvsCrys = getObjectPtr<TH2F>("RawDigitAmpvsCrys");
93  auto RawDigitTimevsCrys = getObjectPtr<TH2F>("RawDigitTimevsCrys");
94  auto hitCrysVsExtrapolatedCrys = getObjectPtr<TH2F>("hitCrysVsExtrapolatedCrys");
95 
99  TH1F* IntegralVsCrysID = new TH1F("IntegralVsCrysID", "Integral of EnVsCrysID for each crystal;crystal ID;Entries", 8736, 0, 8736);
100  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)", 8736, 0,
101  8736);
102  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
103  8736, 0, 8736);
104  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
105  8736, 0, 8736);
106 
107  for (int crysID = 0; crysID < 8736; crysID++) {
108  TH1D* hEnergy = EnVsCrysID->ProjectionY("hEnergy", crysID + 1, crysID + 1);
109  int Integral = hEnergy->Integral();
110  IntegralVsCrysID->SetBinContent(crysID + 1, Integral);
111 
112  double TotEntries = CalibEntriesvsCrys->GetBinContent(crysID + 1);
113 
114  double expectedE = 0.;
115  if (TotEntries > 0.) {expectedE = ExpEvsCrys->GetBinContent(crysID + 1) / TotEntries;}
116  AverageExpECrys->SetBinContent(crysID + 1, expectedE);
117 
118  double calibconst = 0.;
119  if (TotEntries > 0.) {calibconst = ElecCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
120  AverageElecCalib->SetBinContent(crysID + 1, calibconst);
121 
122  calibconst = 0.;
123  if (TotEntries > 0.) {calibconst = InitialCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
124  AverageInitCalib->SetBinContent(crysID + 1, calibconst);
125  }
126 
129  TFile* histfile = new TFile("eclMuMuEAlgorithm.root", "recreate");
130  TrkPerCrysID->Write();
131  EnVsCrysID->Write();
132  IntegralVsCrysID->Write();
133  AverageExpECrys->Write();
134  AverageElecCalib->Write();
135  AverageInitCalib->Write();
136  RawDigitAmpvsCrys->Write();
137  RawDigitTimevsCrys->Write();
138  hitCrysVsExtrapolatedCrys->Write();
139 
142  if (!performFits) {
143  B2RESULT("eclMuMuEAlgorithm has not been asked to perform fits; copying input histograms and quitting");
144  histfile->Close();
145  return c_NotEnoughData;
146  }
147 
150  bool sufficientData = true;
151  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
152  if (IntegralVsCrysID->GetBinContent(crysID + 1) < minEntries) {
153  if (storeConst == 1) {B2RESULT("eclMuMuEAlgorithm: crystal " << crysID << " has insufficient statistics: " << IntegralVsCrysID->GetBinContent(crysID + 1) << ". Requirement is " << minEntries);}
154  sufficientData = false;
155  break;
156  }
157  }
158 
160  if (!sufficientData && storeConst == 1) {
161  histfile->Close();
162  return c_NotEnoughData;
163  }
164 
169  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV", 8736, 0, 8736);
170  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Peak energy (GeV)", 8736, 0, 8736);
171 
173  TH1F* PeakVsCrysID = new TH1F("PeakVsCrysID", "Peak of Novo fit vs crystal ID;crystal ID;Peak normalized energy", 8736, 0,
174  8736);
175  TH1F* effSigVsCrysID = new TH1F("effSigVsCrysID", "effSigma vs crystal ID;crystal ID;sigma)", 8736, 0, 8736);
176  TH1F* etaVsCrysID = new TH1F("etaVsCrysID", "eta vs crystal ID;crystal ID;Novo eta parameter", 8736, 0, 8736);
177  TH1F* constVsCrysID = new TH1F("constVsCrysID", "fit constant vs crystal ID;crystal ID;fit constant", 8736, 0, 8736);
178  TH1F* normVsCrysID = new TH1F("normVsCrysID", "Novosibirsk normalization vs crystal ID;crystal ID;normalization", 8736, 0, 8736);
179  TH1F* fitLimitVsCrysID = new TH1F("fitLimitVsCrysID", "fit range upper limit vs crystal ID;crystal ID;upper fit limit", 8736, 0,
180  8736);
181  TH1F* StatusVsCrysID = new TH1F("StatusVsCrysID", "Fit status vs crystal ID;crystal ID;Fit status", 8736, 0, 8736);
182 
184  TH1F* hStatus = new TH1F("hStatus", "Fit status", 25, -5, 20);
185  TH1F* hPeak = new TH1F("hPeak", "Peaks of normalized energy distributions, successful fits;Peak of Novosibirsk fit", 200, 0.8, 1.2);
186  TH1F* fracPeakUnc = new TH1F("fracPeakUnc", "Fractional uncertainty on peak uncertainty, successful fits;Uncertainty on peak", 100,
187  0, 0.1);
188  TH1F* nIterations = new TH1F("nIterations", "Number of times histogram was fit;Number of iterations", 20, -0.5, 19.5);
189 
190 
193  bool allFitsOK = true;
194  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
195 
197  TString name = "Enormalized";
198  name += crysID;
199  TH1D* hEnergy = EnVsCrysID->ProjectionY(name, crysID + 1, crysID + 1);
200 
202  double histMin = hEnergy->GetXaxis()->GetXmin();
203  double histMax = hEnergy->GetXaxis()->GetXmax();
204  TF1* func = new TF1("eclNovoConst", eclNovoConst, histMin, histMax, 5);
205  func->SetParNames("normalization", "peak", "effSigma", "eta", "const");
206  func->SetParLimits(1, peakMin, peakMax);
207  func->SetParLimits(2, effSigMin, effSigMax);
208  func->SetParLimits(3, etaMin, etaMax);
209  func->SetParLimits(4, constMin, constMax);
210 
212  hEnergy->GetXaxis()->SetRangeUser(peakMin, peakMax);
213  int maxBin = hEnergy->GetMaximumBin();
214  double peakE = hEnergy->GetBinLowEdge(maxBin);
215  double peakEUnc = 0.;
216  double normalization = hEnergy->GetMaximum();
217  double normUnc = 0.;
218  double effSigma = hEnergy->GetRMS();
219  double sigmaUnc = 0.;
220  hEnergy->GetXaxis()->SetRangeUser(histMin, histMax);
221 
223  double fitlow = 0.25;
224  double fithigh = peakE + 2.5 * effSigma;
225 
227  double eta = etaNom;
228  double etaUnc = 0.;
229 
231  int il0 = hEnergy->GetXaxis()->FindBin(fitlow);
232  int il1 = hEnergy->GetXaxis()->FindBin(fitlow + 0.1);
233  double constant = hEnergy->Integral(il0, il1) / (1 + il1 - il0);
234  if (constant < 0.01 * normalization) {constant = 0.01 * normalization;}
235  double constUnc = 0.;
236 
238  double dIter = 0.1 * (histMax - histMin) / hEnergy->GetNbinsX();
239  double fitProb(0.);
240  double highold(0.), higholdold(0.);
241  bool fixConst = false;
242  int nIter = 0;
243  bool fitHist = IntegralVsCrysID->GetBinContent(crysID + 1) >= minEntries; /* fit only if enough events */
244 
247  while (fitHist) {
248 
250  func->SetParameters(normalization, peakE, effSigma, eta, constant);
251  if (fixConst) { func->FixParameter(4, 0); }
252 
254  if (crysID < 1152 || crysID > 7775 || !findExpValues) {func->FixParameter(3, etaNom);}
255 
257  hEnergy->Fit(func, "LIQ", "", fitlow, fithigh);
258  nIter++;
259  fitHist = false;
260  normalization = func->GetParameter(0);
261  normUnc = func->GetParError(0);
262  peakE = func->GetParameter(1);
263  peakEUnc = func->GetParError(1);
264  effSigma = func->GetParameter(2);
265  sigmaUnc = func->GetParError(2);
266  eta = func->GetParameter(3);
267  etaUnc = func->GetParError(3);
268  constant = func->GetParameter(4);
269  constUnc = func->GetParError(4);
270  fitProb = func->GetProb();
271 
273  double peak = func->Eval(peakE) - constant;
274  double tRatio = (func->Eval(fithigh) - constant) / peak;
275  if (tRatio < tRatioMin || tRatio > tRatioMax) {
276  double targetY = constant + 0.5 * (tRatioMin + tRatioMax) * peak;
277  higholdold = highold;
278  highold = fithigh;
279  fithigh = func->GetX(targetY, peakE, histMax);
280  fitHist = true;
281 
283  if (abs(fithigh - higholdold) < dIter) {fithigh = 0.5 * (highold + higholdold); }
284 
286  if (nIter > maxIterations - 3) {fithigh = 0.33333 * (fithigh + highold + higholdold); }
287  }
288 
290  if (constant < constTol && !fixConst) {
291  constant = 0;
292  fixConst = true;
293  fitHist = true;
294  }
295 
297  if (nIter == maxIterations) {fitHist = false;}
298  B2DEBUG(200, crysID << " " << nIter << " " << peakE << " " << constant << " " << tRatio << " " << fithigh);
299  }
300 
303  int iStatus = fitOK; // success
304  if (nIter == maxIterations) {iStatus = iterations;} // too many iterations
305 
307  if (normalization < constRatio * constant) {iStatus = noPeak;}
308 
310  if (fitProb <= minFitLimit || (fitProb < minFitProbIter && iStatus == iterations)) {iStatus = poorFit;}
311 
313  if ((peakE < peakMin + peakTol) || (peakE > peakMax - peakTol)) {iStatus = atLimit;}
314  if ((effSigma < effSigMin + effSigTol) || (effSigma > effSigMax - effSigTol)) {iStatus = atLimit;}
315  if ((eta < etaMin + etaTol) || (eta > etaMax - etaTol)) {iStatus = atLimit;}
316  if (constant > constMax - constTol) {iStatus = atLimit;}
317 
318  //** No fit
319  if (nIter == 0) {iStatus = notFit;} // not fit
320 
322  int histbin = crysID + 1;
323  PeakVsCrysID->SetBinContent(histbin, peakE);
324  PeakVsCrysID->SetBinError(histbin, peakEUnc);
325  effSigVsCrysID->SetBinContent(histbin, effSigma);
326  effSigVsCrysID->SetBinError(histbin, sigmaUnc);
327  etaVsCrysID->SetBinContent(histbin, eta);
328  etaVsCrysID->SetBinError(histbin, etaUnc);
329  constVsCrysID->SetBinContent(histbin, constant);
330  constVsCrysID->SetBinError(histbin, constUnc);
331  normVsCrysID->SetBinContent(histbin, normalization);
332  normVsCrysID->SetBinError(histbin, normUnc);
333  fitLimitVsCrysID->SetBinContent(histbin, fithigh);
334  fitLimitVsCrysID->SetBinError(histbin, 0);
335  StatusVsCrysID->SetBinContent(histbin, iStatus);
336  StatusVsCrysID->SetBinError(histbin, 0);
337 
339  hStatus->Fill(iStatus);
340  nIterations->Fill(nIter);
341  if (iStatus >= iterations) {
342  hPeak->Fill(peakE);
343  fracPeakUnc->Fill(peakEUnc / peakE);
344  }
345 
347  B2INFO("cellID " << crysID + 1 << " status = " << iStatus << " fit probability = " << fitProb);
348  histfile->cd();
349  hEnergy->Write();
350 
351  } /* end of loop over crystals */
352 
355  for (int crysID = 0; crysID < 8736; crysID++) {
356  int histbin = crysID + 1;
357  double fitstatus = StatusVsCrysID->GetBinContent(histbin);
358  double peakE = PeakVsCrysID->GetBinContent(histbin);
359  double fracpeakEUnc = PeakVsCrysID->GetBinError(histbin) / peakE;
360 
362  if (fitstatus < iterations) {
363  peakE = -1.;
364  fracpeakEUnc = 0.;
365  if (histbin >= cellIDLo && histbin <= cellIDHi) {
366  B2RESULT("eclMuMuEAlgorithm: cellID " << histbin << " is not a successful fit. Status = " << fitstatus);
367  allFitsOK = false;
368  }
369  }
370 
372  if (findExpValues) {
373  double inputExpE = abs(AverageExpECrys->GetBinContent(histbin));
374  ExpEnergyperCrys->SetBinContent(histbin, inputExpE * peakE);
375  ExpEnergyperCrys->SetBinError(histbin, fracpeakEUnc * inputExpE * peakE);
376  } else {
377 
379  double inputCalib = abs(AverageInitCalib->GetBinContent(histbin));
380  CalibVsCrysID->SetBinContent(histbin, inputCalib / peakE);
381  CalibVsCrysID->SetBinError(histbin, fracpeakEUnc * inputCalib / peakE);
382  }
383  }
384 
387  bool DBsuccess = false;
388  if (storeConst == 0 || (storeConst == 1 && allFitsOK)) {
389  DBsuccess = true;
390  if (findExpValues) {
391 
393  std::vector<float> tempE;
394  std::vector<float> tempUnc;
395  for (int crysID = 0; crysID < 8736; crysID++) {
396  tempE.push_back(ExpEnergyperCrys->GetBinContent(crysID + 1));
397  tempUnc.push_back(ExpEnergyperCrys->GetBinError(crysID + 1));
398  }
399  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
400  ExpectedE->setCalibVector(tempE, tempUnc);
401  saveCalibration(ExpectedE, "ECLExpMuMuE");
402  B2RESULT("eclCosmicEAlgorithm: successfully stored expected energies ECLExpMuMuE");
403 
404  } else {
405 
407  std::vector<float> tempCalib;
408  std::vector<float> tempCalibUnc;
409  for (int crysID = 0; crysID < 8736; crysID++) {
410  tempCalib.push_back(CalibVsCrysID->GetBinContent(crysID + 1));
411  tempCalibUnc.push_back(CalibVsCrysID->GetBinError(crysID + 1));
412  }
413  ECLCrystalCalib* MuMuECalib = new ECLCrystalCalib();
414  MuMuECalib->setCalibVector(tempCalib, tempCalibUnc);
415  saveCalibration(MuMuECalib, "ECLCrystalEnergyMuMu");
416  B2RESULT("eclMuMuEAlgorithm: successfully stored ECLCrystalEnergyMuMu calibration constants");
417  }
418  }
419 
423  PeakVsCrysID->Write();
424  effSigVsCrysID->Write();
425  etaVsCrysID->Write();
426  constVsCrysID->Write();
427  normVsCrysID->Write();
428  fitLimitVsCrysID->Write();
429  StatusVsCrysID->Write();
430  hPeak->Write();
431  fracPeakUnc->Write();
432  nIterations->Write();
433  hStatus->Write();
434 
436  if (findExpValues) {
437  ExpEnergyperCrys->Write();
438  } else {
439  CalibVsCrysID->Write();
440  }
441  histfile->Close();
442 
445  dummy = (TH1F*)gROOT->FindObject("PeakVsCrysID"); delete dummy;
446  dummy = (TH1F*)gROOT->FindObject("effSigVsCrysID"); delete dummy;
447  dummy = (TH1F*)gROOT->FindObject("etaVsCrysID"); delete dummy;
448  dummy = (TH1F*)gROOT->FindObject("constVsCrysID"); delete dummy;
449  dummy = (TH1F*)gROOT->FindObject("normVsCrysID"); delete dummy;
450  dummy = (TH1F*)gROOT->FindObject("fitLimitVsCrysID"); delete dummy;
451  dummy = (TH1F*)gROOT->FindObject("StatusVsCrysID"); delete dummy;
452  dummy = (TH1F*)gROOT->FindObject("fitProbSame"); delete dummy;
453  dummy = (TH1F*)gROOT->FindObject("fracPeakUnc"); delete dummy;
454  dummy = (TH1F*)gROOT->FindObject("nIterations"); delete dummy;
455  dummy = (TH1F*)gROOT->FindObject("hStatus"); delete dummy;
456  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
457  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
458 
459 
462  if (storeConst == -1) {
463  B2RESULT("eclMuMuEAlgorithm performed fits but was not asked to store contants");
464  return c_Failure;
465  } else if (!DBsuccess) {
466  if (findExpValues) { B2RESULT("eclMuMuEAlgorithm: failed to store expected values"); }
467  else { B2RESULT("eclMuMuEAlgorithm: failed to store calibration constants"); }
468  return c_Failure;
469  }
470  return c_OK;
471 }

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

◆ convertPyExpRun()

ExpRun convertPyExpRun ( PyObject *  pyObj)
inherited

Performs the conversion of PyObject to ExpRun.

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

Definition at line 63 of file CalibrationAlgorithm.cc.

◆ execute()

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

Runs calibration over vector of runs for a given iteration.

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

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

◆ getCollectorName()

std::string getCollectorName ( ) const
inlineinherited

Alias for prefix.

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

Definition at line 174 of file CalibrationAlgorithm.h.

◆ setInputFileNames() [1/2]

void setInputFileNames ( PyObject *  inputFileNames)
inherited

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

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

Definition at line 159 of file CalibrationAlgorithm.cc.

◆ setInputFileNames() [2/2]

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

Set the input file names used for this algorithm.

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

Definition at line 187 of file CalibrationAlgorithm.cc.

Member Data Documentation

◆ cellIDLo

int cellIDLo

..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events

First cellID to be fit

Definition at line 49 of file eclMuMuEAlgorithm.h.

◆ storeConst

int storeConst

controls which values are written to the database.

0 (default): store value found by successful fits, or -|input value| otherwise; -1 : do not store values 1 : store values if every fit for [cellIDLo,cellIDHi] was successful

Definition at line 57 of file eclMuMuEAlgorithm.h.


The documentation for this class was generated from the following files:
prepareAsicCrosstalkSimDB.e
e
aux.
Definition: prepareAsicCrosstalkSimDB.py:53
Belle2::ECL::eclMuMuEAlgorithm::notFit
int notFit
no fit performed
Definition: eclMuMuEAlgorithm.h:73
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::ECL::eclMuMuEAlgorithm::storeConst
int storeConst
controls which values are written to the database.
Definition: eclMuMuEAlgorithm.h:57
Belle2::ECL::eclMuMuEAlgorithm::tRatioMin
double tRatioMin
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclMuMuEAlgorithm.h:53
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::ECLCrystalCalib
General DB object to store one calibration number per ECL crystal.
Definition: ECLCrystalCalib.h:34
Belle2::ECL::eclMuMuEAlgorithm::atLimit
int atLimit
a parameter is at the limit; fit not useable
Definition: eclMuMuEAlgorithm.h:70
Belle2::ECL::eclMuMuEAlgorithm::poorFit
int poorFit
low chi square; fit not useable
Definition: eclMuMuEAlgorithm.h:71
Belle2::ECL::eclMuMuEAlgorithm::findExpValues
bool findExpValues
if true, fits are used to find expected energy deposit for each crystal instead of the calibration co...
Definition: eclMuMuEAlgorithm.h:56
Belle2::ECL::eclMuMuEAlgorithm::cellIDHi
int cellIDHi
Last cellID to be fit.
Definition: eclMuMuEAlgorithm.h:50
Belle2::ECL::eclMuMuEAlgorithm::performFits
bool performFits
if false, input histograms are copied to output, but no fits are done.
Definition: eclMuMuEAlgorithm.h:55
Belle2::ECL::eclMuMuEAlgorithm::cellIDLo
int cellIDLo
..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events
Definition: eclMuMuEAlgorithm.h:49
Belle2::CalibrationAlgorithm::c_Failure
@ c_Failure
Failed =3 in Python.
Definition: CalibrationAlgorithm.h:54
Belle2::ECL::eclMuMuEAlgorithm::noPeak
int noPeak
Novosibirsk component of fit is negligible; fit not useable.
Definition: eclMuMuEAlgorithm.h:72
Belle2::ECL::eclMuMuEAlgorithm::minEntries
int minEntries
All crystals to be fit must have at least minEntries events in the fit range.
Definition: eclMuMuEAlgorithm.h:51
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::ECLCrystalCalib::setCalibVector
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
Definition: ECLCrystalCalib.h:48
Belle2::ECL::eclMuMuEAlgorithm::fitOK
int fitOK
fit is OK
Definition: eclMuMuEAlgorithm.h:68
Belle2::ECL::eclMuMuEAlgorithm::tRatioMax
double tRatioMax
Fit range is adjusted so that fit at upper endpoint is between tRatioMin and tRatioMax of peak.
Definition: eclMuMuEAlgorithm.h:54
Belle2::ECL::eclMuMuEAlgorithm::iterations
int iterations
fit reached max number of iterations, but is useable
Definition: eclMuMuEAlgorithm.h:69
Belle2::ECL::eclMuMuEAlgorithm::maxIterations
int maxIterations
no more than maxIteration iterations
Definition: eclMuMuEAlgorithm.h:52