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
HpsFitResult.h
Go to the documentation of this file.
1#ifndef __HPS_FIT_RESULT_H__
2#define __HPS_FIT_RESULT_H__
3
4#include <TFitResult.h>
5#include <TFitResultPtr.h>
6
9
16
17 public:
18
20 HpsFitResult();
21
23
29 TFitResultPtr getBkgFitResult() { return bkg_result_; };
30
36 TFitResultPtr getBkgToysFitResult() { return bkg_toy_result_; };
37
43 TFitResultPtr getCompFitResult() { return comp_result_; };
44
50 void addLikelihood(double likelihood) { _likelihoods.push_back(likelihood); }
51
57 void addSignalYield(double signal_yield) { _signal_yields.push_back(signal_yield); }
58
64 double getCorrectedMass() const { return _cmass; };
65
71 double getIntegral() { return _integral; };
72
78 double getMass() const { return mass_hypo_; };
79
85 double getQ0() { return q0_; };
86
92 double getPValue() { return p_value_; };
93
99 double getFullBkgRate();
100
106 double getFullBkgRateError();
107
113 float getSignalYield() { return comp_result_->Parameter(poly_order_ + 1); };
114
120 float getSignalYieldErr() { return comp_result_->ParError(poly_order_ + 1); };
121
127 double getUpperLimit() { return upper_limit_; };
128
135
141 double getWindowSize() { return window_size_; };
142
148 std::vector<double> getLikelihoods() { return _likelihoods; }
149
155 std::vector<double> getSignalYields() { return _signal_yields; }
156
157 //-------------//
158 // Setters //
159 //-------------//
160
166 void setIntegral(double integral) { _integral = integral; };
167
173 void setQ0(double q0) { q0_ = q0; };
174
180 void setPValue(double p_value) { p_value_ = p_value; };
181
187 void setBkgFitResult(TFitResultPtr bkg_result) { bkg_result_ = bkg_result; };
188
194 void setBkgToysFitResult(TFitResultPtr bkg_toy_result) { bkg_toy_result_ = bkg_toy_result; };
195
201 void setCompFitResult(TFitResultPtr comp_result) { comp_result_ = comp_result; };
202
208 void setMass(double mass) { mass_hypo_ = mass; };
209
215 void setCorrectedMass(double cmass) { _cmass = cmass; };
216
222 void setPolyOrder(int poly_order) { poly_order_ = poly_order; };
223
229 void setBkgModelType(FitFunction::BkgModel bkg_model) { bkg_model_ = bkg_model; };
230
236 void setUpperLimit(double upper_limit) { upper_limit_ = upper_limit; };
237
243 void setUpperLimitPValue(double upper_limit_p_value) { _upper_limit_p_value = upper_limit_p_value; }
244
250 void setWindowSize(double window_size) { window_size_ = window_size; };
251
257 void setBinWidth(double bin_width) { bin_width_ = bin_width; };
258
259 private:
261 TFitResultPtr bkg_result_{nullptr};
262
264 TFitResultPtr bkg_toy_result_{nullptr};
265
267 TFitResultPtr comp_result_{nullptr};
268
270 std::vector<double> _likelihoods;
271
273 std::vector<double> _signal_yields;
274
276 double _integral{0};
277
279 double mass_hypo_;
280
282 double _cmass{0};
283
285 double q0_;
286
288 double p_value_;
289
291 double poly_order_{0};
292
295
298
300 double _upper_limit_p_value{-9999};
301
302 /*** Size of the fit window. */
304
305 /*** Size of the fit window. */
307
308}; // HpsFitResult
309
310#endif // __HPS_FIT_RESULT_H__
BkgModel
description
Definition FitFunction.h:39
description
void setCompFitResult(TFitResultPtr comp_result)
Set the result from the signal+background only fit.
void setBinWidth(double bin_width)
Set the size of the fit window used to get this results.
double bin_width_
double getIntegral()
Get the integral within the fit window.
void setIntegral(double integral)
Set the integral.
TFitResultPtr getBkgFitResult()
Get background fit result.
void setBkgModelType(FitFunction::BkgModel bkg_model)
Sets the type of background fit function used by the fitter.
void setQ0(double q0)
Set Q0.
TFitResultPtr comp_result_
TFitResultPtr bkg_toy_result_
double poly_order_
TFitResultPtr getCompFitResult()
Get the complete fit result.
void addLikelihood(double likelihood)
description
double window_size_
float getSignalYieldErr()
Get the error of the signal yield.
double mass_hypo_
double getUpperLimit()
Get the upper fit limit.
void setPolyOrder(int poly_order)
Set the order polynomial used by the fitter.
double upper_limit_
void setCorrectedMass(double cmass)
Set the corrected mass.
double getUpperLimitPValue()
Get the p-value at the upper limit.
double _upper_limit_p_value
double getWindowSize()
Get the size of the fit window.
void addSignalYield(double signal_yield)
description
double getCorrectedMass() const
Get the corrected mass.
void setPValue(double p_value)
Set the p-value.
void setMass(double mass)
Set mass hypothesis used for this fit.
void setBkgToysFitResult(TFitResultPtr bkg_toy_result)
Set the result from the background only fit for generating toys.
double getPValue()
description
std::vector< double > _likelihoods
std::vector< double > getLikelihoods()
Get the likelihoods.
double getFullBkgRate()
Get the background rate obtained from the signal+background hit at the mass hypo.
void setUpperLimit(double upper_limit)
Set the 2 sigma upper limit.
TFitResultPtr bkg_result_
std::vector< double > _signal_yields
void setBkgFitResult(TFitResultPtr bkg_result)
Set the result from the background only fit.
FitFunction::BkgModel bkg_model_
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.
void setUpperLimitPValue(double upper_limit_p_value)
Set the p-value at the end of the upper limit calculation.
TFitResultPtr getBkgToysFitResult()
Get background toy model fit result.
std::vector< double > getSignalYields()
Get the signal yields.
void setWindowSize(double window_size)
Set the size of the fit window used to get this results.
double getMass() const
Get the mass hypothesis for the fit.