|
Laura++
3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
|
Class for defining the rescattering model.
More...
#include <LauRescatteringRes.hh>
|
| LauRescatteringRes (LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters) |
| Constructor. More...
|
|
virtual | ~LauRescatteringRes () |
| Destructor.
|
|
virtual void | initialise () |
| Initialise the model.
|
|
virtual LauComplex | amplitude (const LauKinematics *kinematics) |
| Get the complex dynamical amplitude. More...
|
|
virtual LauAbsResonance::LauResonanceModel | getResonanceModel () const |
| Get the resonance model type. More...
|
|
virtual void | setResonanceParameter (const TString &name, const Double_t value) |
| Set value of the various parameters. More...
|
|
virtual void | floatResonanceParameter (const TString &name) |
| Allow the various parameters to float in the fit. More...
|
|
virtual LauParameter * | getResonanceParameter (const TString &name) |
| Access the given resonance parameter. More...
|
|
virtual const std::vector< LauParameter * > & | getFloatingParameters () |
| Retrieve the resonance parameters, e.g. so that they can be loaded into a fit. More...
|
|
| LauAbsResonance (LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters) |
| Constructor (for use by standard resonances) More...
|
|
| LauAbsResonance (const TString &resName, const Int_t resPairAmpInt, const LauDaughters *daughters, const Int_t resSpin) |
| Constructor (for use by K-matrix components) More...
|
|
virtual | ~LauAbsResonance () |
| Destructor.
|
|
LauSpinType | getSpinType () const |
| Get the spin type. More...
|
|
const TString & | getResonanceName () const |
| Get the name of the resonance. More...
|
|
const TString & | getSanitisedName () const |
| Get the name of the resonance. More...
|
|
Int_t | getPairInt () const |
| Get the integer to identify which DP axis the resonance belongs to. More...
|
|
Int_t | getSpin () const |
| Get the spin of the resonance. More...
|
|
Int_t | getCharge () const |
| Get the charge of the resonance. More...
|
|
Double_t | getMass () const |
| Get the mass of the resonance. More...
|
|
Double_t | getWidth () const |
| Get the width of the resonance. More...
|
|
LauParameter * | getMassPar () |
| Get the mass parameter of the resonance. More...
|
|
LauParameter * | getWidthPar () |
| Get the width parameter of the resonance. More...
|
|
virtual Bool_t | preSymmetrised () const |
| Is the amplitude pre-symmetrised? More...
|
|
Bool_t | flipHelicity () const |
| Get the helicity flip flag. More...
|
|
void | flipHelicity (const Bool_t boolean) |
| Set the helicity flip flag. More...
|
|
Bool_t | ignoreMomenta () const |
| Get the ignore momenta flag. More...
|
|
void | ignoreMomenta (const Bool_t boolean) |
| Set the ignore momenta flag. More...
|
|
Bool_t | ignoreSpin () const |
| Get the ignore spin flag. More...
|
|
void | ignoreSpin (const Bool_t boolean) |
| Set the ignore spin flag. More...
|
|
Bool_t | ignoreBarrierScaling () const |
| Get the ignore barrier factor scaling flag. More...
|
|
void | ignoreBarrierScaling (const Bool_t boolean) |
| Set the ignore barrier factor scaling flag. More...
|
|
void | changeResonance (const Double_t newMass, const Double_t newWidth, const Int_t newSpin) |
| Allow the mass, width and spin of the resonance to be changed. More...
|
|
void | changeBWBarrierRadii (const Double_t resRadius, const Double_t parRadius) |
| Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed. More...
|
|
void | fixMass (const Bool_t parFixed) |
| Fix or release the resonance mass. More...
|
|
void | fixWidth (const Bool_t parFixed) |
| Fix or release the resonance width. More...
|
|
Bool_t | fixMass () const |
| Get the status of resonance mass (fixed or released) More...
|
|
Bool_t | fixWidth () const |
| Get the status of resonance width (fixed or released) More...
|
|
void | setSpinType (const LauSpinType spinType) |
| Set the spin formalism to be used. More...
|
|
void | setBarrierRadii (LauBlattWeisskopfFactor *resFactor, LauBlattWeisskopfFactor *parFactor) |
| Set the form factor model and parameters. More...
|
|
void | fixBarrierRadii (const Bool_t fixResRadius, const Bool_t fixParRadius) |
| Fix or release the Blatt-Weisskopf barrier radii.
|
|
Bool_t | fixResRadius () const |
| Get the status of resonance barrier radius (fixed or released)
|
|
Bool_t | fixParRadius () const |
| Get the status of parent barrier radius (fixed or released)
|
|
Double_t | getResRadius () const |
| Get the radius of the resonance barrier factor.
|
|
Double_t | getParRadius () const |
| Get the radius of the parent barrier factor.
|
|
|
enum | LauResonanceModel {
BW,
RelBW,
GS,
Flatte,
Sigma,
Kappa,
Dabba,
LASS,
LASS_BW,
LASS_NR,
EFKLLM,
KMatrix,
FlatNR,
NRModel,
BelleNR,
PowerLawNR,
BelleSymNR,
BelleSymNRNoInter,
TaylorNR,
PolNR,
Pole,
PolarFFNR,
PolarFFSymNR,
PolarFFSymNRNoInter,
Rescattering,
Rescattering2,
RescatteringNoInter,
MIPW_MagPhase,
MIPW_RealImag,
GaussIncoh,
RhoOmegaMix_GS,
RhoOmegaMix_RBW,
RhoOmegaMix_GS_1,
RhoOmegaMix_RBW_1
} |
| Define the allowed resonance types. More...
|
|
enum | LauSpinType { Zemach_P,
Zemach_Pstar,
Covariant,
Legendre
} |
| Define the allowed spin formalisms. More...
|
|
static bool | isIncoherentModel (LauResonanceModel model) |
| Is the resonance model incoherent? More...
|
|
Class for defining the rescattering model.
Defines the Rescatering models from PiPi-KK Inelastic Scatering : 2005: J.R. Pelaez, F. J. Ynduráin: PHYSICAL REVIEW D 71, 074016 (2005) 2015: J. H. Alvarenga Nogueira, I. Bediaga, A. B. R. Cavalcante, T. Frederico, and O. Lourenço: PHYSICAL REVIEW D 92, 054010 (2015) 2018: J.R. Pelaez, A.Rodas: Unpublished yet PiPi -> KK scattering up to 1.47 GeV with hyperbolic dispersion relations.
Definition at line 48 of file LauRescatteringRes.hh.
◆ LauRescatteringRes()
Constructor.
- Parameters
-
[in] | resInfo | the object containing information on the resonance name, mass, width, spin, charge, etc. |
[in] | resType | the model of the resonance |
[in] | resPairAmpInt | the number of the daughter not produced by the resonance |
[in] | daughters | the daughter particles |
Definition at line 41 of file LauRescatteringRes.cc.
◆ amplitude()
Get the complex dynamical amplitude.
- Parameters
-
[in] | kinematics | the kinematic variables of the current event |
- Returns
- the complex amplitude
Reimplemented from LauAbsResonance.
Definition at line 154 of file LauRescatteringRes.cc.
◆ fixC0()
Bool_t LauRescatteringRes::fixC0 |
( |
| ) |
const |
|
inlineprotected |
See if the C0 parameter is fixed or floating.
- Returns
- kTRUE if the C0 parameter is fixed, kFALSE otherwise
Definition at line 252 of file LauRescatteringRes.hh.
◆ fixEps1()
Bool_t LauRescatteringRes::fixEps1 |
( |
| ) |
const |
|
inlineprotected |
See if the Eps1 parameter is fixed or floating.
- Returns
- kTRUE if the Eps1 parameter is fixed, kFALSE otherwise
Definition at line 216 of file LauRescatteringRes.hh.
◆ fixEps2()
Bool_t LauRescatteringRes::fixEps2 |
( |
| ) |
const |
|
inlineprotected |
See if the Eps2 parameter is fixed or floating.
- Returns
- kTRUE if the Eps2 parameter is fixed, kFALSE otherwise
Definition at line 234 of file LauRescatteringRes.hh.
◆ fixLambdaKK()
Bool_t LauRescatteringRes::fixLambdaKK |
( |
| ) |
const |
|
inlineprotected |
See if the lambdaKK parameter is fixed or floating.
- Returns
- kTRUE if the lambdaKK parameter is fixed, kFALSE otherwise
Definition at line 144 of file LauRescatteringRes.hh.
◆ fixLambdaPiPi()
Bool_t LauRescatteringRes::fixLambdaPiPi |
( |
| ) |
const |
|
inlineprotected |
See if the lambdaPiPi parameter is fixed or floating.
- Returns
- kTRUE if the lambdaPiPi parameter is fixed, kFALSE otherwise
Definition at line 126 of file LauRescatteringRes.hh.
◆ fixMf()
Bool_t LauRescatteringRes::fixMf |
( |
| ) |
const |
|
inlineprotected |
See if the Mf parameter is fixed or floating.
- Returns
- kTRUE if the Mf parameter is fixed, kFALSE otherwise
Definition at line 180 of file LauRescatteringRes.hh.
◆ fixMprime()
Bool_t LauRescatteringRes::fixMprime |
( |
| ) |
const |
|
inlineprotected |
See if the Mprime parameter is fixed or floating.
- Returns
- kTRUE if the Mprime parameter is fixed, kFALSE otherwise
Definition at line 198 of file LauRescatteringRes.hh.
◆ fixMs()
Bool_t LauRescatteringRes::fixMs |
( |
| ) |
const |
|
inlineprotected |
See if the Ms parameter is fixed or floating.
- Returns
- kTRUE if the Ms parameter is fixed, kFALSE otherwise
Definition at line 162 of file LauRescatteringRes.hh.
◆ floatResonanceParameter()
void LauRescatteringRes::floatResonanceParameter |
( |
const TString & |
name | ) |
|
|
virtual |
Allow the various parameters to float in the fit.
- Parameters
-
[in] | name | the name of the parameter to be floated |
Reimplemented from LauAbsResonance.
Definition at line 307 of file LauRescatteringRes.cc.
◆ getC0()
Double_t LauRescatteringRes::getC0 |
( |
| ) |
const |
|
inlineprotected |
◆ getEps1()
Double_t LauRescatteringRes::getEps1 |
( |
| ) |
const |
|
inlineprotected |
◆ getEps2()
Double_t LauRescatteringRes::getEps2 |
( |
| ) |
const |
|
inlineprotected |
◆ getFloatingParameters()
const std::vector< LauParameter * > & LauRescatteringRes::getFloatingParameters |
( |
| ) |
|
|
virtual |
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
- Returns
- floating parameters of the resonance
Reimplemented from LauAbsResonance.
Definition at line 234 of file LauRescatteringRes.cc.
◆ getLambdaKK()
Double_t LauRescatteringRes::getLambdaKK |
( |
| ) |
const |
|
inlineprotected |
Get the lambdaKK, the term for the KK.
- Returns
- lambdaKK, the term for the KK
Definition at line 138 of file LauRescatteringRes.hh.
◆ getLambdaPiPi()
Double_t LauRescatteringRes::getLambdaPiPi |
( |
| ) |
const |
|
inlineprotected |
Get the lambdaPiPi, the term for the PiPi.
- Returns
- lambdaPiPi, the term for the PiPi
Definition at line 120 of file LauRescatteringRes.hh.
◆ getMf()
Double_t LauRescatteringRes::getMf |
( |
| ) |
const |
|
inlineprotected |
◆ getMprime()
Double_t LauRescatteringRes::getMprime |
( |
| ) |
const |
|
inlineprotected |
◆ getMs()
Double_t LauRescatteringRes::getMs |
( |
| ) |
const |
|
inlineprotected |
◆ getResonanceModel()
◆ getResonanceParameter()
LauParameter * LauRescatteringRes::getResonanceParameter |
( |
const TString & |
name | ) |
|
|
virtual |
Access the given resonance parameter.
- Parameters
-
[in] | name | the name of the parameter |
- Returns
- the corresponding parameter
Reimplemented from LauAbsResonance.
Definition at line 381 of file LauRescatteringRes.cc.
◆ resAmp()
LauComplex LauRescatteringRes::resAmp |
( |
Double_t |
mass, |
|
|
Double_t |
spinTerm |
|
) |
| |
|
protectedvirtual |
Complex resonant amplitude.
- Parameters
-
[in] | mass | appropriate invariant mass for the resonance |
[in] | spinTerm | Zemach spin term |
Implements LauAbsResonance.
Definition at line 226 of file LauRescatteringRes.cc.
◆ setC0()
void LauRescatteringRes::setC0 |
( |
const Double_t |
C0 | ) |
|
|
protected |
◆ setEps1()
void LauRescatteringRes::setEps1 |
( |
const Double_t |
Eps1 | ) |
|
|
protected |
◆ setEps2()
void LauRescatteringRes::setEps2 |
( |
const Double_t |
Eps2 | ) |
|
|
protected |
◆ setLambdaKK()
void LauRescatteringRes::setLambdaKK |
( |
const Double_t |
lambda | ) |
|
|
protected |
Set the parameter lambdaKK, the term for the KK.
- Parameters
-
[in] | lambda | the term for the KK |
Definition at line 418 of file LauRescatteringRes.cc.
◆ setLambdaPiPi()
void LauRescatteringRes::setLambdaPiPi |
( |
const Double_t |
lambda | ) |
|
|
protected |
Set the parameter lambdaPiPi, the term for the PiPi.
- Parameters
-
[in] | lambda | the term for the PiPi |
Definition at line 411 of file LauRescatteringRes.cc.
◆ setMf()
void LauRescatteringRes::setMf |
( |
const Double_t |
Mf | ) |
|
|
protected |
◆ setMprime()
void LauRescatteringRes::setMprime |
( |
const Double_t |
Mprime | ) |
|
|
protected |
◆ setMs()
void LauRescatteringRes::setMs |
( |
const Double_t |
Ms | ) |
|
|
protected |
◆ setResonanceParameter()
void LauRescatteringRes::setResonanceParameter |
( |
const TString & |
name, |
|
|
const Double_t |
value |
|
) |
| |
|
virtual |
Set value of the various parameters.
- Parameters
-
[in] | name | the name of the parameter to be changed |
[in] | value | the new parameter value |
Reimplemented from LauAbsResonance.
Definition at line 266 of file LauRescatteringRes.cc.
The documentation for this class was generated from the following files:
|