hpstr
The Heavy Photon Search Toolkit for Reconstruction (hpstr) provides an interface to physics data from the HPS experiment saved in the LCIO format and converts it into an ROOT based format.
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
BhMassResSystematicsProcessor.cxx
Go to the documentation of this file.
1#include <math.h>
3
5 : Processor(name, process) {
6 }
7
9
11 std::cout << "Configuring BhMassResSystematicsProcessor" << std::endl;
12 std::string temp_file_name = "";
13 std::string temp_func_name = "";
14 try {
15 debug_ = parameters.getInteger("debug");
16 massSpectrum_ = parameters.getString("massSpectrum");
17 mass_hypo_ = parameters.getDouble("mass_hypo");
18 win_factor_ = parameters.getInteger("win_factor");
19 poly_order_ = parameters.getInteger("poly_order");
20 seed_ = parameters.getInteger("seed");
21 bkg_model_ = parameters.getInteger("bkg_model");
22 res_runs_ = parameters.getInteger("num_iters");
23 res_width_ = parameters.getDouble("res_sigma");
24 temp_file_name = parameters.getString("function_file");
25 function_name_ = parameters.getString("function_name");
26 toy_res_runs_ = parameters.getInteger("toy_num_iters");
27 nToys_ = parameters.getInteger("num_toys");
28 } catch(std::runtime_error& error) {
29 std::cout << error.what() << std::endl;
30 }
31
32 std::cout << "[ DEBUG ] :: Creating ROOT file \"" << temp_file_name << "\"." << std::endl;
33 if(temp_file_name != "") { function_file_ = new TFile(temp_file_name.c_str()); }
34}
35
36void BhMassResSystematicsProcessor::initialize(std::string inFilename, std::string outFilename) {
37 // Initiate the input file.
38 inF_ = new TFile(inFilename.c_str());
39
40 // Get the mass spectrum from file.
41 mass_spec_h = (TH1*) inF_->Get(massSpectrum_.c_str());
42
43 // Get the appropriate background model.
44 FitFunction::BkgModel bkg_fit_model;
45 std::cout << "Background Model ID: " << bkg_model_ << std::endl;
46 switch(bkg_model_) {
47 case 0: bkg_fit_model = FitFunction::BkgModel::EXP_CHEBYSHEV;
48 break;
49 case 1: bkg_fit_model = FitFunction::BkgModel::EXP_CHEBYSHEV;
50 break;
51 case 2: bkg_fit_model = FitFunction::BkgModel::LEGENDRE;
52 break;
53 case 3: bkg_fit_model = FitFunction::BkgModel::EXP_LEGENDRE;
54 break;
55 default: bkg_fit_model = FitFunction::BkgModel::EXP_LEGENDRE;
56 }
57
58 // Initiate the bump hunter.
61 bump_hunter_->setBounds(mass_spec_h->GetXaxis()->GetBinUpEdge(mass_spec_h->FindFirstBinAbove()),
62 mass_spec_h->GetXaxis()->GetBinLowEdge(mass_spec_h->FindLastBinAbove()));
63
64 // Get the mass resolution parameterization, if it is defined.
65 std::cout << "[ DEBUG ] :: Attempting to acquire parameterization function." << std::endl;
66 if(function_file_ != nullptr && function_name_ != "") {
67 std::cout << "[ DEBUG ] :: Function file and name are defined!" << std::endl;
68 std::cout << "[ DEBUG ] :: Function name is \"" << function_name_ << "\"." << std::endl;
69 TF1* massResSysFunc = (TF1*) function_file_->Get(function_name_.c_str());
70 std::cout << "[ DEBUG ] :: Acquired function object." << std::endl;
71 res_width_ = massResSysFunc->Eval(mass_hypo_) / 100.0;
72 std::cout << "[ DEBUG ] :: Resolution width for " << mass_hypo_ << " MeV is " << res_width_ << "." << std::endl;
73 std::cout << "Using mass resolution error parameterization function " << function_name_ << " from file." << std::endl;
74 } else {
75 std::cout << "Using constant mass resolution error " << res_width_ << "." << std::endl;
76 }
77
78 // Enabled debug output if relevant.
79 if(debug_ > 0) { bump_hunter_->enableDebug(); }
80
81 // Initiate the tuple maker.
82 flat_tuple_ = new FlatTupleMaker(outFilename.c_str(), "fit_toys");
83
84 // Setup tuple branches
85 flat_tuple_->addVariable("seed");
86 flat_tuple_->addVariable("poly_order");
87 flat_tuple_->addVariable("win_factor");
88 flat_tuple_->addVariable("bkg_model");
89 flat_tuple_->addVariable("iterations");
90 flat_tuple_->addVariable("res_width");
91 flat_tuple_->addVariable("num_toys");
92 flat_tuple_->addVariable("toy_iterations");
93
94 // Create the nominal value variables.
95 flat_tuple_->addVariable("nominal_bkg_total");
96 flat_tuple_->addVariable("nominal_corr_mass");
97 flat_tuple_->addVariable("nominal_mass_hypo");
98 flat_tuple_->addVariable("nominal_window_size");
99 flat_tuple_->addVariable("nominal_resolution_scale");
100 flat_tuple_->addVariable("nominal_mass_resolution");
101 flat_tuple_->addVariable("nominal_bkg_chi2_prob");
102 flat_tuple_->addVariable("nominal_bkgsig_chi2_prob");
103 flat_tuple_->addVariable("nominal_bkg_edm");
104 flat_tuple_->addVariable("nominal_bkg_minuit_status");
105 flat_tuple_->addVariable("nominal_bkg_nll");
106
107 flat_tuple_->addVariable("nominal_edm");
108 flat_tuple_->addVariable("nominal_minuit_status");
109 flat_tuple_->addVariable("nominal_nll");
110 flat_tuple_->addVariable("nominal_p_value");
111 flat_tuple_->addVariable("nominal_q0");
112 flat_tuple_->addVariable("nominal_bkg_rate_mass_hypo");
113 flat_tuple_->addVariable("nominal_bkg_rate_mass_hypo_err");
114 flat_tuple_->addVariable("nominal_sig_yield");
115 flat_tuple_->addVariable("nominal_sig_yield_err");
116 flat_tuple_->addVariable("nominal_upper_limit");
117 flat_tuple_->addVariable("nominal_ul_p_value");
118
119 // Create the iteration variable vectors.
120 flat_tuple_->addVector("bkg_total");
121 flat_tuple_->addVector("corr_mass");
122 flat_tuple_->addVector("mass_hypo");
123 flat_tuple_->addVector("window_size");
124 flat_tuple_->addVector("resolution_scale");
125 flat_tuple_->addVector("mass_resolution");
126 flat_tuple_->addVector("bkg_chi2_prob");
127 flat_tuple_->addVector("bkgsig_chi2_prob");
128 flat_tuple_->addVector("bkg_edm");
129 flat_tuple_->addVector("bkg_minuit_status");
130 flat_tuple_->addVector("bkg_nll");
131
132 flat_tuple_->addVector("edm");
133 flat_tuple_->addVector("minuit_status");
134 flat_tuple_->addVector("nll");
135 flat_tuple_->addVector("p_value");
136 flat_tuple_->addVector("q0");
137 flat_tuple_->addVector("bkg_rate_mass_hypo");
138 flat_tuple_->addVector("bkg_rate_mass_hypo_err");
139 flat_tuple_->addVector("sig_yield");
140 flat_tuple_->addVector("sig_yield_err");
141 flat_tuple_->addVector("upper_limit");
142 flat_tuple_->addVector("ul_p_value");
143}
144
146 std::cout << "Running on mass spectrum: " << massSpectrum_ << std::endl;
147 std::cout << "Running with polynomial order: " << poly_order_ << std::endl;
148 std::cout << "Running with window factor: " << win_factor_ << std::endl;
149 std::cout << "Running on mass hypothesis: " << mass_hypo_ << std::endl;
150
151 // Set the tuple values which are the same for all runs.
156 flat_tuple_->setVariableValue("iterations", res_runs_);
158 flat_tuple_->setVariableValue("num_toys", nToys_);
159 flat_tuple_->setVariableValue("toy_iterations", toy_res_runs_);
160
161 // Run one fit with the nominal resolution.
162 HpsFitResult* nominal_result = bump_hunter_->performSearch(mass_spec_h, mass_hypo_, false, false);
163 TFitResultPtr nominal_bkg_result = nominal_result->getBkgFitResult();
164 TFitResultPtr nominal_sig_result = nominal_result->getCompFitResult();
165
166 // Write the invariant mass histogram to the output file.
167 mass_spec_h->Write();
168
169 // Set the tuple values.
170 flat_tuple_->setVariableValue("nominal_bkg_total", nominal_result->getIntegral());
171 flat_tuple_->setVariableValue("nominal_corr_mass", nominal_result->getCorrectedMass());
172 flat_tuple_->setVariableValue("nominal_mass_hypo", nominal_result->getMass());
173 flat_tuple_->setVariableValue("nominal_window_size", nominal_result->getWindowSize());
174 flat_tuple_->setVariableValue("nominal_resolution_scale", 1.00);
176
177 flat_tuple_->setVariableValue("nominal_bkg_chi2_prob", nominal_bkg_result->Prob());
178 flat_tuple_->setVariableValue("nominal_bkgsig_chi2_prob", nominal_sig_result->Prob());
179 flat_tuple_->setVariableValue("nominal_bkg_edm", nominal_bkg_result->Edm());
180 flat_tuple_->setVariableValue("nominal_bkg_minuit_status", nominal_bkg_result->Status());
181 flat_tuple_->setVariableValue("nominal_bkg_nll", nominal_bkg_result->MinFcnValue());
182
183 flat_tuple_->setVariableValue("nominal_edm", nominal_sig_result->Edm());
184 flat_tuple_->setVariableValue("nominal_minuit_status", nominal_sig_result->Status());
185 flat_tuple_->setVariableValue("nominal_nll", nominal_sig_result->MinFcnValue());
186 flat_tuple_->setVariableValue("nominal_p_value", nominal_result->getPValue());
187 flat_tuple_->setVariableValue("nominal_q0", nominal_result->getQ0());
188 flat_tuple_->setVariableValue("nominal_bkg_rate_mass_hypo", nominal_result->getFullBkgRate());
189 flat_tuple_->setVariableValue("nominal_bkg_rate_mass_hypo_err", nominal_result->getFullBkgRateError());
190 flat_tuple_->setVariableValue("nominal_sig_yield", nominal_result->getSignalYield());
191 flat_tuple_->setVariableValue("nominal_sig_yield_err", nominal_result->getSignalYieldErr());
192 flat_tuple_->setVariableValue("nominal_upper_limit", nominal_result->getUpperLimit());
193 flat_tuple_->setVariableValue("nominal_ul_p_value", nominal_result->getUpperLimitPValue());
194
195 // Make a random number generator.
196 TRandom3 *rng = new TRandom3();
197
198 // Iterate for a number of loops equal to the desired number of
199 // resolution scale samples.
200 for(int i = 0; i < res_runs_; i++) {
201 // Set the resolution scale to a random value according to a
202 // normal distribution around the actual resolution with a
203 // width as specified by the user.
204 double res_scale = rng->Gaus(1.00, res_width_);
206
207 // Search for a resonance at the given mass hypothesis.
209
210 // Get the result of the background fit.
211 TFitResultPtr bkg_result = result->getBkgFitResult();
212
213 std::cout << "Filling Fit Results " << std::endl;
214 // Get the result of the signal+background fit.
215 TFitResultPtr sig_result = result->getCompFitResult();
216
217 // Set the fit parameters in the tuple.
218 flat_tuple_->addToVector("bkg_total", result->getIntegral());
219 flat_tuple_->addToVector("corr_mass", result->getCorrectedMass());
220 flat_tuple_->addToVector("mass_hypo", result->getMass());
221 flat_tuple_->addToVector("window_size", result->getWindowSize());
222 flat_tuple_->addToVector("resolution_scale", res_scale);
224
225 // Set the fit results in the tuple.
226 flat_tuple_->addToVector("bkg_chi2_prob", bkg_result->Prob());
227 flat_tuple_->addToVector("bkgsig_chi2_prob", sig_result->Prob());
228 flat_tuple_->addToVector("bkg_edm", bkg_result->Edm());
229 flat_tuple_->addToVector("bkg_minuit_status", bkg_result->Status());
230 flat_tuple_->addToVector("bkg_nll", bkg_result->MinFcnValue());
231
232 flat_tuple_->addToVector("edm", sig_result->Edm());
233 flat_tuple_->addToVector("minuit_status", sig_result->Status());
234 flat_tuple_->addToVector("nll", sig_result->MinFcnValue());
235 flat_tuple_->addToVector("p_value", result->getPValue());
236 flat_tuple_->addToVector("q0", result->getQ0());
237 flat_tuple_->addToVector("bkg_rate_mass_hypo", result->getFullBkgRate());
238 flat_tuple_->addToVector("bkg_rate_mass_hypo_err", result->getFullBkgRateError());
239 flat_tuple_->addToVector("sig_yield", result->getSignalYield());
240 flat_tuple_->addToVector("sig_yield_err", result->getSignalYieldErr());
241 flat_tuple_->addToVector("upper_limit", result->getUpperLimit());
242 flat_tuple_->addToVector("ul_p_value", result->getUpperLimitPValue());
243 }
244
245 // Get the quantile vectors.
246 std::vector<double> vec_sigYield = flat_tuple_->getVector("sig_yield");
247 std::vector<double> vec_sigYieldErr = flat_tuple_->getVector("sig_yield_err");
248 std::vector<double> vec_bkgRate = flat_tuple_->getVector("bkg_rate_mass_hypo");
249 std::vector<double> vec_upperLimit = flat_tuple_->getVector("upper_limit");
250 std::vector<double> vec_pValue = flat_tuple_->getVector("p_value");
251
252 // Sort the vectors.
253 std::sort(vec_sigYield.begin(), vec_sigYield.end());
254 std::sort(vec_sigYieldErr.begin(), vec_sigYieldErr.end());
255 std::sort(vec_bkgRate.begin(), vec_bkgRate.end());
256 std::sort(vec_upperLimit.begin(), vec_upperLimit.end());
257 std::sort(vec_pValue.begin(), vec_pValue.end());
258
259 // Calculate the median and 1 sigma indices.
260 int median = (int) rint(vec_sigYield.size() / 2.0);
261 int sigma_1_u = (int) rint((vec_sigYield.size() / 2.0) + (0.341 * vec_sigYield.size()));
262 int sigma_1_l = (int) rint((vec_sigYield.size() / 2.0) - (0.341 * vec_sigYield.size()));
263
264 // Create variables in the tuple to store the quantile values.
265 flat_tuple_->addVariable("sig_yield_median");
266 flat_tuple_->addVariable("sig_yield_sigma1u");
267 flat_tuple_->addVariable("sig_yield_sigma1l");
268
269 flat_tuple_->addVariable("sig_yield_err_median");
270 flat_tuple_->addVariable("sig_yield_err_sigma1u");
271 flat_tuple_->addVariable("sig_yield_err_sigma1l");
272
273 flat_tuple_->addVariable("bkg_rate_mass_hypo_median");
274 flat_tuple_->addVariable("bkg_rate_mass_hypo_sigma1u");
275 flat_tuple_->addVariable("bkg_rate_mass_hypo_sigma1l");
276
277 flat_tuple_->addVariable("upper_limit_median");
278 flat_tuple_->addVariable("upper_limit_sigma1u");
279 flat_tuple_->addVariable("upper_limit_sigma1l");
280
281 flat_tuple_->addVariable("p_value_median");
282 flat_tuple_->addVariable("p_value_sigma1u");
283 flat_tuple_->addVariable("p_value_sigma1l");
284
285 // Store the quantile values.
286 flat_tuple_->setVariableValue("sig_yield_median", vec_sigYield[median]);
287 flat_tuple_->setVariableValue("sig_yield_sigma1u", vec_sigYield[sigma_1_u]);
288 flat_tuple_->setVariableValue("sig_yield_sigma1l", vec_sigYield[sigma_1_l]);
289
290 flat_tuple_->setVariableValue("sig_yield_err_median", vec_sigYieldErr[median]);
291 flat_tuple_->setVariableValue("sig_yield_err_sigma1u", vec_sigYieldErr[sigma_1_u]);
292 flat_tuple_->setVariableValue("sig_yield_err_sigma1l", vec_sigYieldErr[sigma_1_l]);
293
294 flat_tuple_->setVariableValue("bkg_rate_mass_hypo_median", vec_bkgRate[median]);
295 flat_tuple_->setVariableValue("bkg_rate_mass_hypo_sigma1u", vec_bkgRate[sigma_1_u]);
296 flat_tuple_->setVariableValue("bkg_rate_mass_hypo_sigma1l", vec_bkgRate[sigma_1_l]);
297
298 flat_tuple_->setVariableValue("upper_limit_median", vec_upperLimit[median]);
299 flat_tuple_->setVariableValue("upper_limit_sigma1u", vec_upperLimit[sigma_1_u]);
300 flat_tuple_->setVariableValue("upper_limit_sigma1l", vec_upperLimit[sigma_1_l]);
301
302 flat_tuple_->setVariableValue("p_value_median", vec_pValue[median]);
303 flat_tuple_->setVariableValue("p_value_sigma1u", vec_pValue[sigma_1_u]);
304 flat_tuple_->setVariableValue("p_value_sigma1l", vec_pValue[sigma_1_l]);
305
306
307 // Set the parameters for toy generation. Most of these should be
308 // defaults. Systematic studies should not require the more
309 // advanced options.
310 int seed_ = 0;
311 int toy_sig_samples_ = 0;
312 double bkg_mult_ = 1.00;
313 TH1* signal_shape_h_ = nullptr;
314
315 // Generate the vector entries for the toy results. A
316 // different set of tuples must be generated for each
317 // toy distribution.
318 flat_tuple_->addVector("toy_bkg_total");
319 flat_tuple_->addVector("toy_corr_mass");
320 flat_tuple_->addVector("toy_mass_hypo");
321 flat_tuple_->addVector("toy_window_size");
322 flat_tuple_->addVector("toy_resolution_scale");
323 flat_tuple_->addVector("toy_mass_resolution");
324 flat_tuple_->addVector("toy_bkg_chi2_prob");
325 flat_tuple_->addVector("toy_bkgsig_chi2_prob");
326 flat_tuple_->addVector("toy_bkg_edm");
327 flat_tuple_->addVector("toy_bkg_minuit_status");
328 flat_tuple_->addVector("toy_bkg_nll");
329
330 flat_tuple_->addVector("toy_edm");
331 flat_tuple_->addVector("toy_minuit_status");
332 flat_tuple_->addVector("toy_nll");
333 flat_tuple_->addVector("toy_p_value");
334 flat_tuple_->addVector("toy_q0");
335 flat_tuple_->addVector("toy_bkg_rate_mass_hypo");
336 flat_tuple_->addVector("toy_bkg_rate_mass_hypo_err");
337 flat_tuple_->addVector("toy_sig_yield");
338 flat_tuple_->addVector("toy_sig_yield_err");
339 flat_tuple_->addVector("toy_upper_limit");
340 flat_tuple_->addVector("toy_ul_p_value");
341
342 // Generate the toys and perform the toy fits.
343 if(nToys_ > 0) {
344 // Generate toy models.
345 std::cout << "Generating " << nToys_ << " Toys" << std::endl;
346 std::vector<TH1*> toys_hist = bump_hunter_->generateToys(mass_spec_h, nToys_, seed_, toy_sig_samples_, bkg_mult_, signal_shape_h_);
347
348 // Perform the fits for each toy model.
349 int toyFitN = 0;
350 for(TH1* hist : toys_hist) {
351 // Store the results of the mass resolution variance for
352 // this toy.
353 std::cout << "Fitting Toy " << toyFitN << std::endl;
354
355 // Fit the toy multiple times with different mass
356 // resolution values.
357 for(int i = 0; i < toy_res_runs_; i++) {
358 // Set the resolution scale to a random value according to a
359 // normal distribution around the actual resolution with a
360 // width as specified by the user.
361 double res_scale = rng->Gaus(1.00, res_width_);
363
364 // Perform the resonance search.
365 HpsFitResult* result = bump_hunter_->performSearch(hist, mass_hypo_, false, false);
366
367 // Get the result of the background fit.
368 TFitResultPtr bkg_result = result->getBkgFitResult();
369
370 // Get the result of the signal+background fit.
371 TFitResultPtr sig_result = result->getCompFitResult();
372
373 // Set the fit parameters in the tuple.
374 flat_tuple_->addToVector("toy_bkg_total", result->getIntegral());
375 flat_tuple_->addToVector("toy_corr_mass", result->getCorrectedMass());
376 flat_tuple_->addToVector("toy_mass_hypo", result->getMass());
377 flat_tuple_->addToVector("toy_window_size", result->getWindowSize());
378 flat_tuple_->addToVector("toy_resolution_scale", res_scale);
380
381 // Set the fit results in the tuple.
382 flat_tuple_->addToVector("toy_bkg_chi2_prob", bkg_result->Prob());
383 flat_tuple_->addToVector("toy_bkgsig_chi2_prob", sig_result->Prob());
384 flat_tuple_->addToVector("toy_bkg_edm", bkg_result->Edm());
385 flat_tuple_->addToVector("toy_bkg_minuit_status", bkg_result->Status());
386 flat_tuple_->addToVector("toy_bkg_nll", bkg_result->MinFcnValue());
387
388 flat_tuple_->addToVector("toy_edm", sig_result->Edm());
389 flat_tuple_->addToVector("toy_minuit_status", sig_result->Status());
390 flat_tuple_->addToVector("toy_nll", sig_result->MinFcnValue());
391 flat_tuple_->addToVector("toy_p_value", result->getPValue());
392 flat_tuple_->addToVector("toy_q0", result->getQ0());
393 flat_tuple_->addToVector("toy_bkg_rate_mass_hypo", result->getFullBkgRate());
394 flat_tuple_->addToVector("toy_bkg_rate_mass_hypo_err", result->getFullBkgRateError());
395 flat_tuple_->addToVector("toy_sig_yield", result->getSignalYield());
396 flat_tuple_->addToVector("toy_sig_yield_err", result->getSignalYieldErr());
397 flat_tuple_->addToVector("toy_upper_limit", result->getUpperLimit());
398 flat_tuple_->addToVector("toy_ul_p_value", result->getUpperLimitPValue());
399 }
400
401 // Fill the tuple.
402 flat_tuple_->fill();
403
404 // Step to the next toy.
405 toyFitN++;
406 }
407 }
408
409 // Fill and write the tuple
411
412 return true;
413}
414
416 inF_->Close();
417 delete inF_;
418 delete bump_hunter_;
419}
420
#define DECLARE_PROCESSOR(CLASS)
Macro which allows the framework to construct a producer given its name during configuration.
Definition Processor.h:139
Insert description here. more details.
TFile * function_file_
The file that contains the mass resolution error parameterization.
virtual void configure(const ParameterSet &parameters)
Configure using given parameters.
int bkg_model_
What background model type to use.
double res_width_
The width of the resolution Gaussian.
int res_runs_
How many resolution variance runs to make.
double seed_
The toy generator seed. This is always zero.
BhMassResSystematicsProcessor(const std::string &name, Process &process)
Constructor.
double mass_hypo_
The signal hypothesis to use in the fit.
int poly_order_
Order of polynomial used to model the background.
BumpHunter * bump_hunter_
The bump hunter manager.
int nToys_
How many toys to generated.
int win_factor_
The factor that determines the size of the mass window as.
virtual void initialize(std::string inFilename, std::string outFilename)
Initialize processor.
int toy_res_runs_
How many resolution variance runs to make.
FlatTupleMaker * flat_tuple_
The flat tuple manager.
TH1 * mass_spec_h
The mass spectrum to fit.
bool asymptotic_limit_
Whether to use the asymptotic upper limit or the power constrained. Defaults to asymptotic.
std::string function_name_
The name of the function object in the error file.
std::string massSpectrum_
The name of the mass spectrum to fit.
description
Definition BumpHunter.h:49
void setWindowSizeUsesResScale(bool window_use_res_scale)
Sets whether the window size is scaled according to the resolution scaling factor.
Definition BumpHunter.h:177
double getMassResolution(double mass)
Definition BumpHunter.h:202
void setResolutionScale(double res_scale)
Set the resolution after instantiation.
Definition BumpHunter.h:170
void enableDebug(bool debug=true)
Enable/disable debug.
Definition BumpHunter.h:117
void setBounds(double low_bound, double high_bound)
Set the histogram bounds.
std::vector< TH1 * > generateToys(TH1 *histogram, double n_toys, int seed, int toy_sig_samples, int bkg_mult=1, TH1 *signal_hist=nullptr)
description
HpsFitResult * performSearch(TH1 *histogram, double mass_hypothesis, bool skip_bkg_fit, bool skip_ul)
Perform a search for a resonance at the given mass hypothesis.
BkgModel
description
Definition FitFunction.h:39
description
std::vector< double > getVector(std::string variable_name)
description
void addVector(std::string vector_name)
description
void addVariable(std::string variable_name)
description
void close()
description
void addToVector(std::string variable_name, double value)
description
void setVariableValue(std::string variable_name, double value)
description
void fill()
description
description
double getIntegral()
Get the integral within the fit window.
TFitResultPtr getBkgFitResult()
Get background fit result.
TFitResultPtr getCompFitResult()
Get the complete fit result.
float getSignalYieldErr()
Get the error of the signal yield.
double getUpperLimit()
Get the upper fit limit.
double getUpperLimitPValue()
Get the p-value at the upper limit.
double getWindowSize()
Get the size of the fit window.
double getCorrectedMass() const
Get the corrected mass.
double getPValue()
description
double getFullBkgRate()
Get the background rate obtained from the signal+background hit at the mass hypo.
float getSignalYield()
Get the signal yield obtained from the signal+background fit.
double getQ0()
description
double getFullBkgRateError()
Get the background rate error from the signal+background fit at the mass hypo.
double getMass() const
Get the mass hypothesis for the fit.
description
Base class for all event processing components.
Definition Processor.h:34