Belle II Software  release-08-01-10
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 nToRebin
 If fewer entries than this, rebin and fix eta parameter.
 
double tRatioMin
 entries/peak at low edge of fit must be greater than this

 
double tRatioMax
 entries/peak at high edge of fit must be greater than this
 
double lowerEdgeThresh
 Lower edge is where the fit = lowerEdgeThresh * peak value.
 
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
 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

int fitOK = 16
 fit is OK
 
int iterations = 8
 fit reached max number of iterations, but is useable
 
int noLowerEdge = 5
 could not determine lower edge of fit
 
int atLimit = 4
 a parameter is at the limit; fit not useable
 
int poorFit = 3
 low chi square; 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)
 allExpRun
 

Detailed Description

Calibrate ecl crystals using muon pair events.

Definition at line 22 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 40 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


Write out the job parameters


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

eta is nominal values


Fit

Set the initial parameters

Fit

bins on either side of this value

look for the target value between these two points


Fit status

did not find lower edge

poor fit

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 and edge 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 65 of file eclMuMuEAlgorithm.cc.

66 {
69  double limitTol = 0.0005; /*< tolerance for checking if a parameter is at the limit */
70  double minFitLimit = 1e-25; /*< cut off for labeling a fit as poor */
71  double peakMin(0.4), peakMax(2.2); /*< range for peak of normalized energy distribution */
72  double peakTol = limitTol * (peakMax - peakMin); /*< fit is at limit if it is within peakTol of min or max */
73  double effSigMin(0.02), effSigMax(0.2); /*< range for effective sigma of normalized energy distribution */
74  double effSigTol = limitTol * (effSigMax - effSigMin); /*< fit is at limit if it is within effSigTol of min or max */
75  double etaNom(-0.41); /*< Nominal tail parameter; fixed to this value for low statistics fits */
76  double etaMin(-1.), etaMax(0.); /*< Novosibirsk tail parameter range */
77  double etaTol = limitTol * (etaMax - etaMin); /*< fit is at limit if it is within etaTol of min or max */
78 
80  gROOT->SetBatch();
81 
84  B2INFO("eclMuMuEAlgorithm parameters:");
85  B2INFO("cellIDLo = " << cellIDLo);
86  B2INFO("cellIDHi = " << cellIDHi);
87  B2INFO("minEntries = " << minEntries);
88  B2INFO("tRatioMin = " << tRatioMin);
89  B2INFO("tRatioMax = " << tRatioMax);
90  B2INFO("lowerEdgeThresh = " << lowerEdgeThresh);
91  B2INFO("performFits = " << performFits);
92  B2INFO("findExpValues = " << findExpValues);
93  B2INFO("storeConst = " << storeConst);
94 
97  TH1F* dummy;
98  dummy = (TH1F*)gROOT->FindObject("IntegralVsCrysID");
99  if (dummy) {delete dummy;}
100  dummy = (TH1F*)gROOT->FindObject("AverageExpECrys");
101  if (dummy) {delete dummy;}
102  dummy = (TH1F*)gROOT->FindObject("AverageElecCalib");
103  if (dummy) {delete dummy;}
104  dummy = (TH1F*)gROOT->FindObject("AverageInitCalib");
105  if (dummy) {delete dummy;}
106 
109  auto TrkPerCrysID = getObjectPtr<TH1F>("TrkPerCrysID");
110  auto EnVsCrysID = getObjectPtr<TH2F>("EnVsCrysID");
111  auto ExpEvsCrys = getObjectPtr<TH1F>("ExpEvsCrys");
112  auto ElecCalibvsCrys = getObjectPtr<TH1F>("ElecCalibvsCrys");
113  auto InitialCalibvsCrys = getObjectPtr<TH1F>("InitialCalibvsCrys");
114  auto CalibEntriesvsCrys = getObjectPtr<TH1F>("CalibEntriesvsCrys");
115  auto RawDigitAmpvsCrys = getObjectPtr<TH2F>("RawDigitAmpvsCrys");
116  auto RawDigitTimevsCrys = getObjectPtr<TH2F>("RawDigitTimevsCrys");
117  auto hitCrysVsExtrapolatedCrys = getObjectPtr<TH2F>("hitCrysVsExtrapolatedCrys");
118 
122  TH1F* IntegralVsCrysID = new TH1F("IntegralVsCrysID", "Integral of EnVsCrysID for each crystal;crystal ID;Entries",
124  TH1F* AverageExpECrys = new TH1F("AverageExpECrys", "Average expected E per crys from collector;Crystal ID;Energy (GeV)",
127  TH1F* AverageElecCalib = new TH1F("AverageElecCalib", "Average electronics calib const vs crystal;Crystal ID;Calibration constant",
129  TH1F* AverageInitCalib = new TH1F("AverageInitCalib", "Average initial calib const vs crystal;Crystal ID;Calibration constant",
131 
132  for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
133  TH1D* hEnergy = EnVsCrysID->ProjectionY("hEnergy", crysID + 1, crysID + 1);
134  int Integral = hEnergy->Integral();
135  IntegralVsCrysID->SetBinContent(crysID + 1, Integral);
136 
137  double TotEntries = CalibEntriesvsCrys->GetBinContent(crysID + 1);
138 
139  double expectedE = 0.;
140  if (TotEntries > 0.) {expectedE = ExpEvsCrys->GetBinContent(crysID + 1) / TotEntries;}
141  AverageExpECrys->SetBinContent(crysID + 1, expectedE);
142 
143  double calibconst = 0.;
144  if (TotEntries > 0.) {calibconst = ElecCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
145  AverageElecCalib->SetBinContent(crysID + 1, calibconst);
146 
147  calibconst = 0.;
148  if (TotEntries > 0.) {calibconst = InitialCalibvsCrys->GetBinContent(crysID + 1) / TotEntries;}
149  AverageInitCalib->SetBinContent(crysID + 1, calibconst);
150  }
151 
154  TFile* histfile = new TFile("eclMuMuEAlgorithm.root", "recreate");
155  TrkPerCrysID->Write();
156  EnVsCrysID->Write();
157  IntegralVsCrysID->Write();
158  AverageExpECrys->Write();
159  AverageElecCalib->Write();
160  AverageInitCalib->Write();
161  RawDigitAmpvsCrys->Write();
162  RawDigitTimevsCrys->Write();
163  hitCrysVsExtrapolatedCrys->Write();
164 
167  if (!performFits) {
168  B2RESULT("eclMuMuEAlgorithm has not been asked to perform fits; copying input histograms and quitting");
169  histfile->Close();
170  return c_NotEnoughData;
171  }
172 
175  bool sufficientData = true;
176  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
177  if (IntegralVsCrysID->GetBinContent(crysID + 1) < minEntries) {
178  if (storeConst == 1) {B2RESULT("eclMuMuEAlgorithm: crystal " << crysID << " has insufficient statistics: " << IntegralVsCrysID->GetBinContent(crysID + 1) << ". Requirement is " << minEntries);}
179  sufficientData = false;
180  break;
181  }
182  }
183 
185  if (!sufficientData && storeConst == 1) {
186  histfile->Close();
187  return c_NotEnoughData;
188  }
189 
194  TH1F* CalibVsCrysID = new TH1F("CalibVsCrysID", "Calibration constant vs crystal ID;crystal ID;counts per GeV",
196  TH1F* ExpEnergyperCrys = new TH1F("ExpEnergyperCrys", "Expected energy per crystal;Crystal ID;Peak energy (GeV)",
198 
200  TH1F* PeakVsCrysID = new TH1F("PeakVsCrysID", "Peak of Novo fit vs crystal ID;crystal ID;Peak normalized energy",
203  TH1F* EdgeVsCrysID = new TH1F("EdgeVsCrysID", "Lower edge of Novo fit vs crystal ID;crystal ID", ECLElementNumbers::c_NCrystals, 0,
205  TH1F* effSigVsCrysID = new TH1F("effSigVsCrysID", "effSigma vs crystal ID;crystal ID;sigma)", ECLElementNumbers::c_NCrystals, 0,
207  TH1F* etaVsCrysID = new TH1F("etaVsCrysID", "eta vs crystal ID;crystal ID;Novo eta parameter", ECLElementNumbers::c_NCrystals, 0,
209  TH1F* normVsCrysID = new TH1F("normVsCrysID", "Novosibirsk normalization vs crystal ID;crystal ID;normalization",
211  TH1F* lowerLimitVsCrysID = new TH1F("lowerLimitVsCrysID", "fit range lower limit vs crystal ID;crystal ID;lower fit limit",
214  TH1F* fitLimitVsCrysID = new TH1F("fitLimitVsCrysID", "fit range upper limit vs crystal ID;crystal ID;upper fit limit",
217  TH1F* StatusVsCrysID = new TH1F("StatusVsCrysID", "Fit status vs crystal ID;crystal ID;Fit status", ECLElementNumbers::c_NCrystals,
219 
221  TH1F* hStatus = new TH1F("hStatus", "Fit status", 25, -5, 20);
222  TH1F* hPeak = new TH1F("hPeak", "Peaks of normalized energy distributions, successful fits;Peak of Novosibirsk fit", 200, 0.8, 1.2);
223  TH1F* fracPeakUnc = new TH1F("fracPeakUnc", "Fractional uncertainty on peak uncertainty, successful fits;Uncertainty on peak", 100,
224  0, 0.1);
225 
226 
229  bool allFitsOK = true;
230  for (int crysID = cellIDLo - 1; crysID < cellIDHi; crysID++) {
231 
233  TString name = "Enormalized";
234  name += crysID;
235  TH1D* hEnergy = EnVsCrysID->ProjectionY(name, crysID + 1, crysID + 1);
236 
238  double histMin = hEnergy->GetXaxis()->GetXmin();
239  double histMax = hEnergy->GetXaxis()->GetXmax();
240  TF1* func = new TF1("eclNovoConst", eclNovoConst, histMin, histMax, 5);
241  func->SetParNames("normalization", "peak", "effSigma", "eta", "const");
242  func->SetParLimits(1, peakMin, peakMax);
243  func->SetParLimits(2, effSigMin, effSigMax);
244  func->SetParLimits(3, etaMin, etaMax);
245 
246  //..Currently not using the constant term
247  double constant = 0.;
248  func->FixParameter(4, constant);
249 
250  //..If low statistics, rebin, and fix eta
251  if (hEnergy->GetEntries() < nToRebin) {
252  hEnergy->Rebin(2);
253  func->FixParameter(3, etaNom);
254  }
255 
257  hEnergy->GetXaxis()->SetRangeUser(peakMin, peakMax);
258  double peak = hEnergy->GetMaximum();
259  int maxBin = hEnergy->GetMaximumBin();
260  double peakE = hEnergy->GetBinLowEdge(maxBin);
261  double peakEUnc = 0.;
262  double normalization = hEnergy->GetMaximum();
263  double normUnc = 0.;
264  double effSigma = hEnergy->GetRMS();
265  double sigmaUnc = 0.;
266  hEnergy->GetXaxis()->SetRangeUser(histMin, histMax);
267  double fitProb = 0.;
268 
270  double eta = etaNom;
271  double etaUnc = 0.;
272 
273  //..Will find the lower edge of normalized energy at the end of the fit
274  double lowerEnEdge = 0.;
275 
276  //..Fit range from set of bins with sufficient entries.
277  double targetY = tRatioMin * peak;
278  int iLow = maxBin;
279  while (hEnergy->GetBinContent(iLow) > targetY) {iLow--;}
280  double fitlow = hEnergy->GetBinLowEdge(iLow);
281 
282  targetY = tRatioMax * peak;
283  int iHigh = maxBin;
284  while (hEnergy->GetBinContent(iHigh) > targetY) {iHigh++;}
285  double fithigh = hEnergy->GetBinLowEdge(iHigh + 1);
286 
289  bool fitHist = IntegralVsCrysID->GetBinContent(crysID + 1) >= minEntries; /* fit only if enough events */
290  if (fitHist) {
291 
293  func->SetParameters(normalization, peakE, effSigma, eta, constant);
294 
296  hEnergy->Fit(func, "LIQ", "", fitlow, fithigh);
297  normalization = func->GetParameter(0);
298  normUnc = func->GetParError(0);
299  peakE = func->GetParameter(1);
300  peakEUnc = func->GetParError(1);
301  effSigma = func->GetParameter(2);
302  sigmaUnc = func->GetParError(2);
303  eta = func->GetParameter(3);
304  etaUnc = func->GetParError(3);
305  fitProb = func->GetProb();
306 
307  //..Lower edge of the fit function is used to find the calibration constant.
308  // Can now use the peak of the fit, instead of the bin content.
309  peak = func->Eval(peakE);
310  targetY = lowerEdgeThresh * peak;
311 
313  iHigh = hEnergy->GetXaxis()->FindBin(peakE) + 1;
314  iLow = hEnergy->GetXaxis()->FindBin(fitlow);
315  int iLast = iHigh;
316  for (int ibin = iHigh; ibin > iLow; ibin--) {
317  double xc = hEnergy->GetBinCenter(ibin);
318  if (func->Eval(xc) > targetY) {iLast = ibin;}
319  }
320  double xHigh = hEnergy->GetBinCenter(iLast);
321  double xLow = hEnergy->GetBinCenter(iLast - 1);
322 
324  if (func->Eval(xLow) < targetY and func->Eval(xHigh) > targetY) {
325  func->SetNpx(1000);
326  lowerEnEdge = func->GetX(targetY, xLow, xHigh);
327  }
328  }
329 
332  int iStatus = fitOK; // success
333 
335  if (lowerEnEdge < 0.01) {iStatus = noLowerEdge;}
336 
338  if (fitProb <= minFitLimit) {iStatus = poorFit;}
339 
341  if ((peakE < peakMin + peakTol) || (peakE > peakMax - peakTol)) {iStatus = atLimit;}
342  if ((effSigma < effSigMin + effSigTol) || (effSigma > effSigMax - effSigTol)) {iStatus = atLimit;}
343  if ((eta < etaMin + etaTol) || (eta > etaMax - etaTol)) {iStatus = atLimit;}
344 
345  //** No fit
346  if (!fitHist) {iStatus = notFit;} // not fit
347 
349  int histbin = crysID + 1;
350  PeakVsCrysID->SetBinContent(histbin, peakE);
351  PeakVsCrysID->SetBinError(histbin, peakEUnc);
352  EdgeVsCrysID->SetBinContent(histbin, lowerEnEdge);
353  EdgeVsCrysID->SetBinError(histbin, peakEUnc); // approximate
354  effSigVsCrysID->SetBinContent(histbin, effSigma);
355  effSigVsCrysID->SetBinError(histbin, sigmaUnc);
356  etaVsCrysID->SetBinContent(histbin, eta);
357  etaVsCrysID->SetBinError(histbin, etaUnc);
358  normVsCrysID->SetBinContent(histbin, normalization);
359  normVsCrysID->SetBinError(histbin, normUnc);
360  lowerLimitVsCrysID->SetBinContent(histbin, fitlow);
361  lowerLimitVsCrysID->SetBinError(histbin, 0);
362  fitLimitVsCrysID->SetBinContent(histbin, fithigh);
363  fitLimitVsCrysID->SetBinError(histbin, 0);
364  StatusVsCrysID->SetBinContent(histbin, iStatus);
365  StatusVsCrysID->SetBinError(histbin, 0);
366 
368  hStatus->Fill(iStatus);
369  if (iStatus >= iterations) {
370  hPeak->Fill(peakE);
371  fracPeakUnc->Fill(peakEUnc / peakE);
372  }
373 
375  B2INFO("cellID " << crysID + 1 << " status = " << iStatus << " fit probability = " << fitProb);
376  histfile->cd();
377  hEnergy->Write();
378 
379  } /* end of loop over crystals */
380 
383  for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
384  int histbin = crysID + 1;
385  double fitstatus = StatusVsCrysID->GetBinContent(histbin);
386  double peakE = PeakVsCrysID->GetBinContent(histbin);
387  double edge = EdgeVsCrysID->GetBinContent(histbin);
388  double fracpeakEUnc = PeakVsCrysID->GetBinError(histbin) / peakE;
389 
391  if (fitstatus < iterations) {
392  peakE = -1.;
393  edge = -1.;
394  fracpeakEUnc = 0.;
395  if (histbin >= cellIDLo && histbin <= cellIDHi) {
396  B2RESULT("eclMuMuEAlgorithm: cellID " << histbin << " is not a successful fit. Status = " << fitstatus);
397  allFitsOK = false;
398  }
399  }
400 
402  if (findExpValues) {
403  double inputExpE = abs(AverageExpECrys->GetBinContent(histbin));
404  ExpEnergyperCrys->SetBinContent(histbin, inputExpE * edge);
405  ExpEnergyperCrys->SetBinError(histbin, fracpeakEUnc * inputExpE * peakE);
406  } else {
407 
409  double inputCalib = abs(AverageInitCalib->GetBinContent(histbin));
410  CalibVsCrysID->SetBinContent(histbin, inputCalib / edge);
411  CalibVsCrysID->SetBinError(histbin, fracpeakEUnc * inputCalib / peakE);
412  }
413  }
414 
417  bool DBsuccess = false;
418  if (storeConst == 0 || (storeConst == 1 && allFitsOK)) {
419  DBsuccess = true;
420  if (findExpValues) {
421 
423  std::vector<float> tempE;
424  std::vector<float> tempUnc;
425  for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
426  tempE.push_back(ExpEnergyperCrys->GetBinContent(crysID + 1));
427  tempUnc.push_back(ExpEnergyperCrys->GetBinError(crysID + 1));
428  }
429  ECLCrystalCalib* ExpectedE = new ECLCrystalCalib();
430  ExpectedE->setCalibVector(tempE, tempUnc);
431  saveCalibration(ExpectedE, "ECLExpMuMuE");
432  B2RESULT("eclCosmicEAlgorithm: successfully stored expected energies ECLExpMuMuE");
433 
434  } else {
435 
437  std::vector<float> tempCalib;
438  std::vector<float> tempCalibUnc;
439  for (int crysID = 0; crysID < ECLElementNumbers::c_NCrystals; crysID++) {
440  tempCalib.push_back(CalibVsCrysID->GetBinContent(crysID + 1));
441  tempCalibUnc.push_back(CalibVsCrysID->GetBinError(crysID + 1));
442  }
443  ECLCrystalCalib* MuMuECalib = new ECLCrystalCalib();
444  MuMuECalib->setCalibVector(tempCalib, tempCalibUnc);
445  saveCalibration(MuMuECalib, "ECLCrystalEnergyMuMu");
446  B2RESULT("eclMuMuEAlgorithm: successfully stored ECLCrystalEnergyMuMu calibration constants");
447  }
448  }
449 
453  PeakVsCrysID->Write();
454  EdgeVsCrysID->Write();
455  effSigVsCrysID->Write();
456  etaVsCrysID->Write();
457  normVsCrysID->Write();
458  lowerLimitVsCrysID->Write();
459  fitLimitVsCrysID->Write();
460  StatusVsCrysID->Write();
461  hPeak->Write();
462  fracPeakUnc->Write();
463  hStatus->Write();
464 
466  if (findExpValues) {
467  ExpEnergyperCrys->Write();
468  } else {
469  CalibVsCrysID->Write();
470  }
471  histfile->Close();
472 
475  dummy = (TH1F*)gROOT->FindObject("PeakVsCrysID"); delete dummy;
476  dummy = (TH1F*)gROOT->FindObject("EdgeVsCrysID"); delete dummy;
477  dummy = (TH1F*)gROOT->FindObject("effSigVsCrysID"); delete dummy;
478  dummy = (TH1F*)gROOT->FindObject("etaVsCrysID"); delete dummy;
479  dummy = (TH1F*)gROOT->FindObject("normVsCrysID"); delete dummy;
480  dummy = (TH1F*)gROOT->FindObject("lowerLimitVsCrysID"); delete dummy;
481  dummy = (TH1F*)gROOT->FindObject("fitLimitVsCrysID"); delete dummy;
482  dummy = (TH1F*)gROOT->FindObject("StatusVsCrysID"); delete dummy;
483  dummy = (TH1F*)gROOT->FindObject("fitProbSame"); delete dummy;
484  dummy = (TH1F*)gROOT->FindObject("fracPeakUnc"); delete dummy;
485  dummy = (TH1F*)gROOT->FindObject("hStatus"); delete dummy;
486  dummy = (TH1F*)gROOT->FindObject("ExpEnergyperCrys"); delete dummy;
487  dummy = (TH1F*)gROOT->FindObject("CalibVsCrysID"); delete dummy;
488 
489 
492  if (storeConst == -1) {
493  B2RESULT("eclMuMuEAlgorithm performed fits but was not asked to store contants");
494  return c_Failure;
495  } else if (!DBsuccess) {
496  if (findExpValues) { B2RESULT("eclMuMuEAlgorithm: failed to store expected values"); }
497  else { B2RESULT("eclMuMuEAlgorithm: failed to store calibration constants"); }
498  return c_Failure;
499  }
500  return c_OK;
501 }
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
@ c_Failure
Failed =3 in Python.
General DB object to store one calibration number per ECL crystal.
void setCalibVector(const std::vector< float > &CalibConst, const std::vector< float > &CalibConstUnc)
Set vector of constants with uncertainties.
double tRatioMin
entries/peak at low edge of fit must be greater than this
int poorFit
low chi square; fit not useable
int iterations
fit reached max number of iterations, but is useable
int noLowerEdge
could not determine lower edge of fit
int storeConst
controls which values are written to the database.
bool performFits
if false, input histograms are copied to output, but no fits are done.
int cellIDHi
Last cellID to be fit.
double lowerEdgeThresh
Lower edge is where the fit = lowerEdgeThresh * peak value.
int cellIDLo
..Parameters to control Novosibirsk fit to energy deposited in each crystal by mu+mu- events
bool findExpValues
if true, fits are used to find expected energy deposit for each crystal instead of the calibration co...
int minEntries
All crystals to be fit must have at least minEntries events in the fit range.
int nToRebin
If fewer entries than this, rebin and fix eta parameter.
int atLimit
a parameter is at the limit; fit not useable
double tRatioMax
entries/peak at high edge of fit must be greater than this
const int c_NCrystals
Number of crystals.

◆ checkPyExpRun()

bool checkPyExpRun ( PyObject *  pyObj)
inherited

Checks that a PyObject can be successfully converted to an ExpRun type.

Checks if the PyObject can be converted to ExpRun.

Definition at line 28 of file CalibrationAlgorithm.cc.

◆ convertPyExpRun()

ExpRun convertPyExpRun ( PyObject *  pyObj)
inherited

Performs the conversion of PyObject to ExpRun.

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

Definition at line 70 of file CalibrationAlgorithm.cc.

◆ execute()

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.

◆ getCollectorName()

std::string getCollectorName ( ) const
inlineinherited

Alias for prefix.

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

Definition at line 164 of file CalibrationAlgorithm.h.

◆ setInputFileNames() [1/2]

void setInputFileNames ( PyObject *  inputFileNames)
inherited

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

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

Definition at line 166 of file CalibrationAlgorithm.cc.

◆ setInputFileNames() [2/2]

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

Set the input file names used for this algorithm.

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

Definition at line 194 of file CalibrationAlgorithm.cc.

Member Data Documentation

◆ 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 32 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 41 of file eclMuMuEAlgorithm.h.


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