MINT2
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Model_independent Class Reference

#include <Model_independent.h>

Inheritance diagram for Model_independent:
BW_BW ILineshape ILineshape MINT::FitParDependent MINT::IFitParDependent MINT::IFitParDependent MINT::IFitParRegister MINT::PolymorphVector< FitParRef > MINT::IFitParDependent

Public Member Functions

 Model_independent (const AssociatedDecayTree &tree, const std::string &namePrefix="")
 
virtual std::string name () const
 
virtual ~Model_independent ()
 
virtual std::complex< double > getVal (IDalitzEvent &evt)
 
std::complex< double > getInterpolatedVal (const double s_inGeV2)
 
void plot ()
 
- Public Member Functions inherited from BW_BW
virtual double prSq () const
 
virtual double prSqForGofM () const
 
virtual double pABSq ()
 
virtual double prSqMax () const
 
 BW_BW (const AssociatedDecayTree &decay, const std::string &lineshapePrefix="", MINT::MinuitParameterSet *mps=0)
 
 BW_BW (const BW_BW &other)
 
virtual ~BW_BW ()
 
virtual void print (IDalitzEvent &evt, std::ostream &out=std::cout)
 
virtual void print (std::ostream &out=std::cout) const
 
virtual DalitzCoordinate getDalitzCoordinate (double nSigma=3) const
 
virtual MINT::counted_ptr< IGenFctgeneratingFunction () const
 
virtual void setGenerationLimits (double mi, double ma)
 
virtual int twoLPlusOne () const
 
virtual std::complex< double > EvtGenValue (IDalitzEvent &evt)
 
- Public Member Functions inherited from ILineshape
virtual ~ILineshape ()
 
- Public Member Functions inherited from MINT::FitParDependent
virtual unsigned int size () const
 
virtual const FitParRefoperator[] (unsigned int i) const
 
virtual FitParRefoperator[] (unsigned int i)
 
virtual bool changedSinceLastCall () const
 
virtual void rememberFitParValues ()
 
virtual bool registerFitParDependence (const IFitParDependent &fpd)
 
bool registerFitParDependence (const FitParRef &fpr)
 
void removeAllFitParDependencies ()
 
 FitParDependent (IFitParRegister *daddy=0)
 
 FitParDependent (const FitParDependent &other, IFitParRegister *newDaddy=0)
 
void listFitParDependencies (std::ostream &os=std::cout) const
 
- Public Member Functions inherited from MINT::PolymorphVector< FitParRef >
 PolymorphVector ()
 
 PolymorphVector (unsigned int N)
 
 PolymorphVector (unsigned int N, const FitParRef &c)
 
 PolymorphVector (const PolymorphVector &other)
 
 PolymorphVector (const typename std::vector< FitParRef > &other)
 
virtual ~PolymorphVector ()
 
std::vector< FitParRef > & theVector ()
 
const std::vector< FitParRef > & theVector () const
 
FitParRefoperator[] (unsigned int i)
 
const FitParRefoperator[] (unsigned int i) const
 
FitParRefat (unsigned int i)
 
const FitParRefat (unsigned int i) const
 
std::vector< FitParRef >::iterator begin ()
 
std::vector< FitParRef >::const_iterator begin () const
 
std::vector< FitParRef >::iterator end ()
 
std::vector< FitParRef >::const_iterator end () const
 
std::vector< FitParRef >::iterator find (const FitParRef &c)
 
std::vector< FitParRef >::const_iterator find (const FitParRef &c) const
 
FitParReffront ()
 
const FitParReffront () const
 
FitParRefback ()
 
const FitParRefback () const
 
unsigned int size () const
 
bool empty () const
 
void push_back (const FitParRef &c)
 
void pop_back ()
 
void erase (typename std::vector< FitParRef >::iterator pos)
 
void erase (typename std::vector< FitParRef >::iterator first, typename std::vector< FitParRef >::iterator last)
 
PolymorphVector< FitParRef > & operator= (const PolymorphVector< FitParRef > &other)
 
void clear ()
 
void resize (unsigned int N)
 
void resize (unsigned int N, const FitParRef &c)
 
 operator const typename std::vector< FitParRef > & () const
 
 operator typename std::vector< FitParRef > & ()
 
bool operator== (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator!= (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator< (const MINT::PolymorphVector< FitParRef > &v2) const
 
bool operator> (const MINT::PolymorphVector< FitParRef > &v2) const
 

Protected Member Functions

std::vector< double > getBinCenterValues_Re ()
 
std::vector< double > getBinCenterValues_Im ()
 
double Re_Bin (int i) const
 
double Im_Bin (int i) const
 
double Bin_1_Re () const
 
double Bin_1_Im () const
 
double Bin_2_Re () const
 
double Bin_2_Im () const
 
double Bin_3_Re () const
 
double Bin_3_Im () const
 
double Bin_4_Re () const
 
double Bin_4_Im () const
 
double Bin_5_Re () const
 
double Bin_5_Im () const
 
double Bin_6_Re () const
 
double Bin_6_Im () const
 
double Bin_7_Re () const
 
double Bin_7_Im () const
 
double Bin_8_Re () const
 
double Bin_8_Im () const
 
double Bin_9_Re () const
 
double Bin_9_Im () const
 
double Bin_10_Re () const
 
double Bin_10_Im () const
 
- Protected Member Functions inherited from BW_BW
const MINT::MinuitParameterSetgetMinuitParameterSet () const
 
MINT::MinuitParameterSetgetMinuitParameterSet ()
 
const std::string & prefix () const
 
bool substitutePDGForReco () const
 
virtual double twoBody_dgtPsq_in_MumsFrame (double mumsMass, double mA, double mB) const
 
virtual double twoBody_dgtPsq_in_MumsPDGFrame () const
 
virtual double twoBody_recodgtPsq_in_MumsPDGFrame () const
 
virtual double twoBody_dgtPsq_in_MumsRecoFrame ()
 
bool setEventPtr (IDalitzEvent &evt) const
 
IDalitzEventgetEvent () const
 
bool setAllFitParameters ()
 
ResonancePropertiesListresonancePropertiesList () const
 
const ResonancePropertiesresonanceProperties () const
 
virtual const ParticlePropertiesmumsProperties () const
 
virtual ResonancePropertiesFitRefmumsFittableProperties () const
 
virtual int mumsPID () const
 
virtual double mumsMass () const
 
virtual double mumsWidth () const
 
virtual double mumsRadius () const
 
virtual double globalRadius () const
 
virtual double Radius () const
 
double GetAlpha () const
 
virtual int lowestPossibleTwoLPlusOne () const
 
virtual int maxDaughterSpinSum () const
 
virtual int minDaughterSpinSum () const
 
virtual int minDaughterSpinSum2 () const
 
virtual int minDaughterSpinSum3 () const
 
virtual int minDaughterSpinSum4 () const
 
virtual int maxDaughterPairSpinSum (int i, int j) const
 
virtual int minDaughterPairSpinSum (int i, int j) const
 
virtual std::string mumsSpin () const
 
virtual int mumsSpinValue () const
 
virtual int mumsParity () const
 
virtual bool startOfDecayChain () const
 
virtual bool parityConservingL (int L) const
 
virtual double mumsRecoMass2 () const
 
virtual double mumsRecoMass () const
 
virtual MultiQuarkContent mumsQuarkContent () const
 
bool nonResonant () const
 
virtual int numDaughters () const
 
virtual TLorentzVector daughterP4 (int i) const
 
virtual double daughterPDGMass (const int &i) const
 
virtual double daughterWidth (int i) const
 
virtual double daughterRecoMass2 (int i) const
 
virtual double daughterRecoMass (int i) const
 
virtual std::string daughterSpin (int i) const
 
virtual int daughterSpinValue (int i) const
 
virtual MultiQuarkContent daughterQuarkContent (int i) const
 
virtual int daughterP (int i) const
 
virtual int dgtrsInternalParity () const
 
virtual MultiQuarkContent dgtrsQuarkContent () const
 
virtual bool isWeakDecay () const
 
virtual double Fr_PDG_BL ()
 
virtual double Fr_BELLE (double prSquared)
 
virtual double Fr_BELLE_Max ()
 
virtual double Fr ()
 
virtual double FrForGofM ()
 
virtual double FrMax ()
 
virtual double GofM ()
 
virtual std::complex< double > BreitWigner ()
 
virtual void resetInternals ()
 
virtual void resetPDG ()
 

Protected Attributes

int _nBins
 
NamedParameter< double > _binCenters
 
std::vector< double > _binCenterVector
 
NamedParameter< std::string > _interpolationTypeString
 
ROOT::Math::Interpolation::Type _interpolationType
 
ROOT::Math::Interpolator * _interpolator_Re
 
ROOT::Math::Interpolator * _interpolator_Im
 
- Protected Attributes inherited from BW_BW
MINT::MinuitParameterSet_mps
 
std::string _prefix
 
MINT::NamedParameter< int > _normBF
 
bool _useGlobalRadius
 
const AssociatedDecayTree_theDecay
 
double _gen_s_mi
 
double _gen_s_ma
 
ResonancePropertiesList_RPL
 
ResonancePropertiesFitRef_fittableResonancePropertiesPtr
 
MINT::FitParRef_fittableGlobalRadiusPtr
 
- Protected Attributes inherited from MINT::PolymorphVector< FitParRef >
std::vector< FitParRef_vec
 

Detailed Description

Definition at line 21 of file Model_independent.h.

Constructor & Destructor Documentation

◆ Model_independent()

Model_independent::Model_independent ( const AssociatedDecayTree tree,
const std::string &  namePrefix = "" 
)

Definition at line 20 of file Model_independent.cpp.

20  :
21  BW_BW(tree, namePrefix)
22  , _binCenters("BinCenters_inGeV2")
23  , _interpolationTypeString("MI::InterpolationType",(std::string)"kCSPLINE")
24  , _interpolationType(ROOT::Math::Interpolation::kCSPLINE)
25 {
27  _nBins = _binCenterVector.size();
28  if((string)_interpolationTypeString == "kLINEAR")_interpolationType = ROOT::Math::Interpolation::kLINEAR;
29  if((string)_interpolationTypeString == "kPOLYNOMIAL")_interpolationType = ROOT::Math::Interpolation::kPOLYNOMIAL;
30  if((string)_interpolationTypeString == "kAKIMA")_interpolationType = ROOT::Math::Interpolation::kAKIMA;
31  if((string)_interpolationTypeString == "kCSPLINE_PERIODIC")_interpolationType = ROOT::Math::Interpolation::kCSPLINE_PERIODIC;
32  if((string)_interpolationTypeString == "kAKIMA_PERIODIC")_interpolationType = ROOT::Math::Interpolation::kAKIMA_PERIODIC;
33 
34  vector<double> binCenterValues_Re = getBinCenterValues_Re();
35  vector<double> binCenterValues_Im = getBinCenterValues_Im();
36  _interpolator_Re = new ROOT::Math::Interpolator(_binCenterVector, binCenterValues_Re, _interpolationType);
37  _interpolator_Im = new ROOT::Math::Interpolator(_binCenterVector, binCenterValues_Im, _interpolationType);
38 }
std::vector< double > getBinCenterValues_Im()
NamedParameter< std::string > _interpolationTypeString
ROOT::Math::Interpolator * _interpolator_Re
BW_BW(const AssociatedDecayTree &decay, const std::string &lineshapePrefix="", MINT::MinuitParameterSet *mps=0)
Definition: BW_BW.cpp:26
std::vector< double > getBinCenterValues_Re()
std::vector< double > _binCenterVector
NamedParameter< double > _binCenters
const std::vector< T > & getVector() const
ROOT::Math::Interpolator * _interpolator_Im
ROOT::Math::Interpolation::Type _interpolationType

◆ ~Model_independent()

virtual Model_independent::~Model_independent ( )
inlinevirtual

Definition at line 30 of file Model_independent.h.

30  {
31  plot();
32  delete _interpolator_Re;
33  delete _interpolator_Im;
34  }
ROOT::Math::Interpolator * _interpolator_Re
ROOT::Math::Interpolator * _interpolator_Im

Member Function Documentation

◆ Bin_10_Im()

double Model_independent::Bin_10_Im ( ) const
inlineprotected

Definition at line 139 of file Model_independent.h.

139  {
141  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_10_Re()

double Model_independent::Bin_10_Re ( ) const
inlineprotected

Definition at line 136 of file Model_independent.h.

136  {
138  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_1_Im()

double Model_independent::Bin_1_Im ( ) const
inlineprotected

Definition at line 85 of file Model_independent.h.

85  {
87  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_1_Re()

double Model_independent::Bin_1_Re ( ) const
inlineprotected

Definition at line 82 of file Model_independent.h.

82  {
84  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_2_Im()

double Model_independent::Bin_2_Im ( ) const
inlineprotected

Definition at line 91 of file Model_independent.h.

91  {
93  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_2_Re()

double Model_independent::Bin_2_Re ( ) const
inlineprotected

Definition at line 88 of file Model_independent.h.

88  {
90  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_3_Im()

double Model_independent::Bin_3_Im ( ) const
inlineprotected

Definition at line 97 of file Model_independent.h.

97  {
99  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_3_Re()

double Model_independent::Bin_3_Re ( ) const
inlineprotected

Definition at line 94 of file Model_independent.h.

94  {
96  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_4_Im()

double Model_independent::Bin_4_Im ( ) const
inlineprotected

Definition at line 103 of file Model_independent.h.

103  {
105  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_4_Re()

double Model_independent::Bin_4_Re ( ) const
inlineprotected

Definition at line 100 of file Model_independent.h.

100  {
102  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_5_Im()

double Model_independent::Bin_5_Im ( ) const
inlineprotected

Definition at line 109 of file Model_independent.h.

109  {
111  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_5_Re()

double Model_independent::Bin_5_Re ( ) const
inlineprotected

Definition at line 106 of file Model_independent.h.

106  {
108  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_6_Im()

double Model_independent::Bin_6_Im ( ) const
inlineprotected

Definition at line 115 of file Model_independent.h.

115  {
117  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_6_Re()

double Model_independent::Bin_6_Re ( ) const
inlineprotected

Definition at line 112 of file Model_independent.h.

112  {
114  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_7_Im()

double Model_independent::Bin_7_Im ( ) const
inlineprotected

Definition at line 121 of file Model_independent.h.

121  {
123  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_7_Re()

double Model_independent::Bin_7_Re ( ) const
inlineprotected

Definition at line 118 of file Model_independent.h.

118  {
120  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_8_Im()

double Model_independent::Bin_8_Im ( ) const
inlineprotected

Definition at line 127 of file Model_independent.h.

127  {
129  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_8_Re()

double Model_independent::Bin_8_Re ( ) const
inlineprotected

Definition at line 124 of file Model_independent.h.

124  {
126  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_9_Im()

double Model_independent::Bin_9_Im ( ) const
inlineprotected

Definition at line 133 of file Model_independent.h.

133  {
135  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ Bin_9_Re()

double Model_independent::Bin_9_Re ( ) const
inlineprotected

Definition at line 130 of file Model_independent.h.

130  {
132  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ getBinCenterValues_Im()

std::vector< double > Model_independent::getBinCenterValues_Im ( )
protected

Definition at line 48 of file Model_independent.cpp.

48  {
49  vector<double> vec;
50  for (int i= 0; i<_nBins; i++) {
51  vec.push_back(Im_Bin(i));
52  }
53  return vec;
54 }
double Im_Bin(int i) const

◆ getBinCenterValues_Re()

std::vector< double > Model_independent::getBinCenterValues_Re ( )
protected

Definition at line 40 of file Model_independent.cpp.

40  {
41  vector<double> vec;
42  for (int i= 0; i<_nBins; i++) {
43  vec.push_back(Re_Bin(i));
44  }
45  return vec;
46 }
double Re_Bin(int i) const

◆ getInterpolatedVal()

std::complex< double > Model_independent::getInterpolatedVal ( const double  s_inGeV2)

Definition at line 68 of file Model_independent.cpp.

68  {
69 
70  if(_binCenterVector[0]>s_inGeV2)return 0;
71  if(_binCenterVector[_nBins-1]<s_inGeV2)return 0;
72 
73  vector<double> binCenterValues_Re = getBinCenterValues_Re();
74  vector<double> binCenterValues_Im = getBinCenterValues_Im();
75 
77  delete _interpolator_Re;
78  delete _interpolator_Im;
79  _interpolator_Re = new ROOT::Math::Interpolator(_binCenterVector, binCenterValues_Re, _interpolationType);
80  _interpolator_Im = new ROOT::Math::Interpolator(_binCenterVector, binCenterValues_Im, _interpolationType);
81  }
82  return complex<double>(_interpolator_Re->Eval(s_inGeV2),_interpolator_Im->Eval(s_inGeV2));
83 }
virtual bool changedSinceLastCall() const
std::vector< double > getBinCenterValues_Im()
ROOT::Math::Interpolator * _interpolator_Re
std::vector< double > getBinCenterValues_Re()
std::vector< double > _binCenterVector
ROOT::Math::Interpolator * _interpolator_Im
ROOT::Math::Interpolation::Type _interpolationType

◆ getVal()

std::complex< double > Model_independent::getVal ( IDalitzEvent evt)
virtual

Reimplemented from BW_BW.

Definition at line 56 of file Model_independent.cpp.

56  {
57  setEventPtr(evt);
59 
60  double formFactor= 1.;
61  if( _normBF == 1 ) formFactor = Fr();
62  else if( _normBF == 0 ) formFactor = Fr_PDG_BL();
63  else if(_normBF == 2 ) formFactor = Fr_BELLE(0.);
64 
65  return formFactor*getInterpolatedVal(mumsRecoMass2()/(GeV*GeV));
66 }
virtual void resetInternals()
Definition: BW_BW.cpp:1114
virtual double Fr()
Definition: BW_BW.cpp:710
bool setEventPtr(IDalitzEvent &evt) const
Definition: BW_BW.cpp:141
std::complex< double > getInterpolatedVal(const double s_inGeV2)
virtual double Fr_PDG_BL()
Definition: BW_BW.cpp:830
virtual double mumsRecoMass2() const
Definition: BW_BW.cpp:548
MINT::NamedParameter< int > _normBF
Definition: BW_BW.h:66
static const double GeV
virtual double Fr_BELLE(double prSquared)
Definition: BW_BW.cpp:720

◆ Im_Bin()

double Model_independent::Im_Bin ( int  i) const
inlineprotected

Definition at line 68 of file Model_independent.h.

68  {
79  else return 0.;
80  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

◆ name()

virtual std::string Model_independent::name ( ) const
inlinevirtual

Reimplemented from BW_BW.

Definition at line 26 of file Model_independent.h.

26  {
27  return "Model_independent("+_theDecay.oneLiner() +")";
28  }
const AssociatedDecayTree & _theDecay
Definition: BW_BW.h:69
void oneLiner(std::stringstream &seam, int generation=0) const
Definition: DDTree.h:375

◆ plot()

void Model_independent::plot ( )

Definition at line 86 of file Model_independent.cpp.

86  {
87 
88  TCanvas* c = new TCanvas();
89 
90  int n = _nBins;
91  vector<double> binCenterValues_Re = getBinCenterValues_Re();
92  vector<double> binCenterValues_Im = getBinCenterValues_Im();
93 
94  std::vector<double> Re(n), Im(n), mag(n), phase(n), delta(n), mass2(n);
95 
96  for (int i=0; i<n; i++) {
97  Re[i] = binCenterValues_Re[i] ;
98  Im[i] = binCenterValues_Im[i] ;
99  mag[i] = pow(binCenterValues_Re[i],2) + pow(binCenterValues_Im[i],2);
100  delta[i] = 1./2. * atan( (binCenterValues_Im[i] - 1./2.)/binCenterValues_Re[i] ) + TMath::Pi() / 4.;
101  phase[i] = arg(complex<double>(Re[i],Im[i]) );
102  mass2[i] = _binCenterVector[i];
103  }
104 
105  double m2_start = _binCenterVector[0];
106  double step = 0.01;
107  int nSteps = ( (_binCenterVector[n-1]) - m2_start)/step;
108  std::vector<double> Re_itp(nSteps), Im_itp(nSteps), mag_itp(nSteps), phase_itp(nSteps), delta_itp(nSteps), mass2_itp(nSteps);
109 
110  for(int i=0; i<nSteps; i++){
111  mass2_itp[i] = m2_start + i * step;
112  Re_itp[i] = getInterpolatedVal(mass2_itp[i]).real();
113  Im_itp[i] = getInterpolatedVal(mass2_itp[i]).imag() ;
114  mag_itp[i] = pow(Re_itp[i],2) + pow(Im_itp[i],2);
115  phase_itp[i] = arg(complex<double>(Re_itp[i],Im_itp[i]) );
116  delta_itp[i] = 1./2. * atan( (Im_itp[i] - 1./2.)/Re_itp[i] ) + TMath::Pi() / 4.;
117  }
118 
119  TGraph* g = new TGraph(n,&Re[0],&Im[0]);
120  g->SetTitle("Argand; Re(A) ; Im(A)");
121  TGraph* g_itp = new TGraph(nSteps,&Re_itp[0],&Im_itp[0]);
122  g_itp->SetLineColor(kRed);
123  g->Draw("A*");
124  g_itp->Draw("Csame");
125  c->Print(("argand_"+anythingToString(mumsPID())+".eps").c_str());
126  delete g;
127  delete g_itp;
128 
129  g = new TGraph(n,&mass2[0],&mag[0]);
130  g->SetTitle(" ; s[GeV^{2}] ; |A|^{2}");
131  g_itp = new TGraph(nSteps,&mass2_itp[0],&mag_itp[0]);
132  g_itp->SetLineColor(kRed);
133  g->Draw("A*");
134  g_itp->Draw("Csame");
135  c->Print(("mag_"+anythingToString(mumsPID())+".eps").c_str());
136  delete g;
137  delete g_itp;
138 
139  g = new TGraph(n,&mass2[0],&phase[0]);
140  g->SetTitle(" ; s[GeV^{2}] ; arg(A)");
141  g_itp = new TGraph(nSteps,&mass2_itp[0],&phase_itp[0]);
142  g_itp->SetLineColor(kRed);
143  g->Draw("A*");
144  g_itp->Draw("Csame");
145  c->Print(("phase_"+anythingToString(mumsPID())+".eps").c_str());
146  delete g;
147  delete g_itp;
148 
149  g = new TGraph(n,&mass2[0],&delta[0]);
150  g->SetTitle(" ; s[GeV^{2}] ; #delta(A)");
151  g_itp = new TGraph(nSteps,&mass2_itp[0],&delta_itp[0]);
152  g_itp->SetLineColor(kRed);
153  g->Draw("A*");
154  g_itp->Draw("Csame");
155  c->Print(("delta_"+anythingToString(mumsPID())+".eps").c_str());
156  delete g;
157  delete g_itp;
158 }
std::vector< double > getBinCenterValues_Im()
virtual int mumsPID() const
Definition: BW_BW.cpp:517
std::complex< double > getInterpolatedVal(const double s_inGeV2)
std::vector< double > getBinCenterValues_Re()
std::vector< double > _binCenterVector
std::string anythingToString(const T &anything)
Definition: Utils.h:62
static const double g

◆ Re_Bin()

double Model_independent::Re_Bin ( int  i) const
inlineprotected

Definition at line 54 of file Model_independent.h.

54  {
65  else return 0.;
66  }
virtual ResonancePropertiesFitRef & mumsFittableProperties() const
Definition: BW_BW.cpp:414

Member Data Documentation

◆ _binCenters

NamedParameter<double> Model_independent::_binCenters
protected

Definition at line 42 of file Model_independent.h.

◆ _binCenterVector

std::vector<double> Model_independent::_binCenterVector
protected

Definition at line 43 of file Model_independent.h.

◆ _interpolationType

ROOT::Math::Interpolation::Type Model_independent::_interpolationType
protected

Definition at line 49 of file Model_independent.h.

◆ _interpolationTypeString

NamedParameter<std::string> Model_independent::_interpolationTypeString
protected

Definition at line 48 of file Model_independent.h.

◆ _interpolator_Im

ROOT::Math::Interpolator* Model_independent::_interpolator_Im
protected

Definition at line 52 of file Model_independent.h.

◆ _interpolator_Re

ROOT::Math::Interpolator* Model_independent::_interpolator_Re
protected

Definition at line 51 of file Model_independent.h.

◆ _nBins

int Model_independent::_nBins
protected

Definition at line 41 of file Model_independent.h.


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