KaliVeda
Toolkit for HIC analysis
KVEvaporationSpectrum.h
1 #pragma once
52 #include "TF1.h"
53 #include "TMultiGraph.h"
54 
55 class KVEvaporationSpectrum : public TF1
56 {
58  {
59  bool con_sig;
60  double norm_const;
61  spectrum_calculator(bool sig_constrain = true)
62  : con_sig{sig_constrain}
63  {
64  if(con_sig)
66  else
68  }
69  double operator()(double*,double*) const;
70  bool is_sigma_constrained() const { return con_sig; }
71  };
73  {
74  int Nc;
76  multi_comp_spectrum_calculator(int ncomp, bool sig_constrain = true)
77  : Nc{ncomp}, SC{sig_constrain}
78  {}
79  double operator()(double*,double*) const;
80  bool is_sigma_constrained() const { return SC.is_sigma_constrained(); }
81  };
82 
83 public:
85  KVEvaporationSpectrum(int ncomponents=1, bool constrain_sigma = true)
86  : KVEvaporationSpectrum("KVEvaporationSpectrum", ncomponents, constrain_sigma)
87  {}
88  KVEvaporationSpectrum(const TString& name, int ncomponents=1, bool constrain_sigma = true);
89  void SetInitialParameters(const TH1*);
90  auto GetBarrier(int comp = 0) const
91  {
92  return GetParameter(get_npar()*comp+1);
93  }
94  auto GetBarrierUncertainty(int comp = 0) const
95  {
96  return GetParError(get_npar()*comp+1);
97  }
98  auto GetTemperature(int comp = 0) const
99  {
100  return GetParameter(get_npar()*comp+2);
101  }
102  auto GetTemperatureUncertainty(int comp = 0) const
103  {
104  return GetParError(get_npar()*comp+2);
105  }
106 
108  {
109  return GetNpar()/get_npar();
110  }
111  TMultiGraph* GetComponents(const std::vector<Color_t>& line_colors={kGray+3, kGray+2, kGray+1, kGray}) const;
112 
114 
115  private:
116  int get_npar() const
117  {
118  return spe_calc.is_sigma_constrained() ? 3 : 4;
119  }
120 };
121 
int Int_t
kGray
#define ClassDefOverride(name, id)
Fittable parametrisation of energy spectra for evaporated particles.
void SetInitialParameters(const TH1 *)
auto GetTemperature(int comp=0) const
auto GetTemperatureUncertainty(int comp=0) const
multi_comp_spectrum_calculator spe_calc
auto GetBarrier(int comp=0) const
TMultiGraph * GetComponents(const std::vector< Color_t > &line_colors={kGray+3, kGray+2, kGray+1, kGray}) const
auto GetBarrierUncertainty(int comp=0) const
KVEvaporationSpectrum(int ncomponents=1, bool constrain_sigma=true)
virtual Double_t GetParameter(const TString &name) const
virtual Double_t GetParError(Int_t ipar) const
virtual Int_t GetNpar() const
constexpr Double_t E()
Double_t Sqrt(Double_t x)
constexpr Double_t TwoPi()
multi_comp_spectrum_calculator(int ncomp, bool sig_constrain=true)