laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauPolarFormFactorSymNR.hh
Go to the documentation of this file.
1 
2 /*
3 Copyright 2018 University of Warwick
4 
5 Licensed under the Apache License, Version 2.0 (the "License");
6 you may not use this file except in compliance with the License.
7 You may obtain a copy of the License at
8 
9  http://www.apache.org/licenses/LICENSE-2.0
10 
11 Unless required by applicable law or agreed to in writing, software
12 distributed under the License is distributed on an "AS IS" BASIS,
13 WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14 See the License for the specific language governing permissions and
15 limitations under the License.
16 */
17 
18 /*
19 Laura++ package authors:
20 John Back
21 Paul Harrison
22 Thomas Latham
23 */
24 
38 #ifndef LAU_POLAR_FORM_FACTOR_SYM_NR
39 #define LAU_POLAR_FORM_FACTOR_SYM_NR
40 
41 #include "LauAbsResonance.hh"
42 #include "LauComplex.hh"
43 
44 #include "TString.h"
45 
46 class LauKinematics;
47 class LauParameter;
48 
50 
51  public:
53 
61  const Int_t resPairAmpInt,
62  const LauDaughters* daughters );
63 
65  virtual ~LauPolarFormFactorSymNR();
66 
68  virtual void initialise();
69 
71 
75  virtual LauComplex amplitude( const LauKinematics* kinematics );
76 
78 
82 
84 
88  virtual void setResonanceParameter( const TString& name, const Double_t value );
89 
91 
94  virtual void floatResonanceParameter( const TString& name );
95 
97 
101  virtual LauParameter* getResonanceParameter( const TString& name );
102 
104 
107  virtual const std::vector<LauParameter*>& getFloatingParameters();
108 
109  protected:
111 
114  void setLambda( const Double_t lambda );
115 
117 
120  Double_t getLambda() const { return ( lambda_ != 0 ) ? lambda_->value() : 0.0; }
121 
123 
126  Bool_t fixLambda() const { return ( lambda_ != 0 ) ? lambda_->fixed() : kTRUE; }
127 
129 
133  virtual LauComplex resAmp( Double_t mass, Double_t spinTerm );
134 
135  private:
138 
141 
144 
147 
148  ClassDef( LauPolarFormFactorSymNR, 0 )
149 };
150 
151 #endif
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
virtual void initialise()
Initialise the model.
void setLambda(const Double_t lambda)
Set the parameter lambda, the NR shape parameter.
Double_t value() const
The value of the parameter.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
LauPolarFormFactorSymNR & operator=(const LauPolarFormFactorSymNR &rhs)
Copy assignment operator (not implemented)
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
virtual ~LauPolarFormFactorSymNR()
Destructor.
Class for defining a nonresonant form factor model.
File containing declaration of LauAbsResonance class.
Double_t getLambda() const
Get the parameter lambda, the NR shape parameter.
virtual LauAbsResonance::LauResonanceModel getResonanceModel() const
Get the resonance model type.
Bool_t fixLambda() const
See if the lambda parameter is fixed or floating.
Class for defining a complex number.
Definition: LauComplex.hh:61
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
File containing declaration of LauComplex class.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauPolarFormFactorSymNR(LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
LauParameter * lambda_
The NR shape parameter.
Class for defining the properties of a resonant particle.
LauAbsResonance::LauResonanceModel model_
The model to use.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Get the complex dynamical amplitude.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
Class for calculating 3-body kinematic quantities.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
LauPolarFormFactorSymNR(const LauPolarFormFactorSymNR &rhs)
Copy constructor (not implemented)
LauResonanceModel
Define the allowed resonance types.