laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauPolarFormFactorSymNR.cc
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 
30 
31 #include "LauConstants.hh"
32 #include "LauDaughters.hh"
33 #include "LauParameter.hh"
34 #include "LauResonanceInfo.hh"
35 
36 #include "TMath.h"
37 
38 #include <iostream>
39 
42  const Int_t resPairAmpInt,
43  const LauDaughters* daughters ) :
44  LauAbsResonance( resInfo, resPairAmpInt, daughters ),
45  lambda_( 0 ),
46  model_( resType )
47 
48 {
49  TString parName = this->getSanitisedName();
50  parName += "_lambda";
51  lambda_ = resInfo->getExtraParameter( parName );
52  if ( lambda_ == 0 ) {
53  lambda_ = new LauParameter( parName, 1.0, 0.0, 10.0, kTRUE );
54  lambda_->secondStage( kTRUE );
55  resInfo->addExtraParameter( lambda_ );
56  }
57 }
58 
60 {
61 }
62 
64 {
65  const LauDaughters* daughters = this->getDaughters();
66  if ( ! daughters->gotSymmetricalDP() ) {
67  std::cerr << "WARNING in LauPolarFormFactorSymNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate."
68  << std::endl;
69  }
70 
71  Int_t resPairAmpInt = this->getPairInt();
72  if ( resPairAmpInt == 3 ) {
73  std::cerr << "WARNING in LauPolarFormFactorSymNR::initialise : This lineshape is intended to be on the symmetrised axes of the DP."
74  << std::endl;
75  }
76 
79  std::cerr << "WARNING in LauPolarFormFactorSymNR::initialise : Unknown model requested, defaulting to Polar Form Factor."
80  << std::endl;
82  }
83 
84  if ( ( model_ != LauAbsResonance::PolarFFSymNR ) && ( this->getSpin() != 0 ) ) {
85  std::cerr << "WARNING in LauPolarFormFactorSymNR::initialise : Non-zero spin will be ignored for this model - perhaps you should use LauAbsResonance::PolarFFSymNRNoInter instead"
86  << std::endl;
87  }
88 
89  // NB we do not need to call setSpinType(LauAbsResonance::Legendre) here (as is done in LauPolarFormFactorNR) since override the amplitude method and explicitly use calcLegendrePoly
90 }
91 
92 LauComplex LauPolarFormFactorSymNR::resAmp( Double_t mass, Double_t spinTerm )
93 {
94  std::cerr << "ERROR in LauPolarFormFactorSymNR : This method should never be called."
95  << std::endl;
96  std::cerr << " : Returning zero amplitude for mass = " << mass
97  << " and spinTerm = " << spinTerm << "." << std::endl;
98  return LauComplex( 0.0, 0.0 );
99 }
100 
102 {
103  // This function returns the complex dynamical amplitude for a Polar Form Factor Non-Resonant distribution
104 
105  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
106  // Have s<->t symmetry already done in Dynamics flip function.
107  // For Kpipi or similar plots, one can use the separate terms
108  // and consider them as two separate components with their own mag and phase.
109  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
110  // individual amplitudes (with the same value of lambda).
111 
112  // Calculate Mandelstam variables.
113  // s = m_13^2, t = m_23^2
114  const Double_t s = kinematics->getm13Sq();
115  const Double_t t = kinematics->getm23Sq();
116 
117  Double_t magnitude( 1.0 );
118 
119  const Double_t lambda = this->getLambda();
120 
122 
123  magnitude = 1.0 / ( 1.0 + s / ( lambda * lambda ) ) + 1.0 / ( 1.0 + t / ( lambda * lambda ) );
124 
126 
127  magnitude = ( s <= t ) ? 1.0 / ( 1.0 + s / ( lambda * lambda ) )
128  : 1.0 / ( 1.0 + t / ( lambda * lambda ) );
129  }
130 
131  LauComplex resAmplitude( magnitude, 0.0 );
132 
133  return resAmplitude;
134 }
135 
136 const std::vector<LauParameter*>& LauPolarFormFactorSymNR::getFloatingParameters()
137 {
138  this->clearFloatingParameters();
139 
140  if ( ! this->fixLambda() ) {
141  this->addFloatingParameter( lambda_ );
142  }
143 
144  return this->getParameters();
145 }
146 
147 void LauPolarFormFactorSymNR::setResonanceParameter( const TString& name, const Double_t value )
148 {
149  // Set various parameters for the lineshape
150  if ( name == "lambda" ) {
151  this->setLambda( value );
152  std::cout << "INFO in LauPolarFormFactorSymNR::setResonanceParameter : Setting parameter lambda = "
153  << this->getLambda() << std::endl;
154  } else {
155  std::cerr << "WARNING in LauPolarFormFactorSymNR::setResonanceParameter: Parameter name not reconised. No parameter changes made."
156  << std::endl;
157  }
158 }
159 
161 {
162  if ( name == "lambda" ) {
163  if ( lambda_->fixed() ) {
164  lambda_->fixed( kFALSE );
165  this->addFloatingParameter( lambda_ );
166  } else {
167  std::cerr << "WARNING in LauPolarFormFactorSymNR::floatResonanceParameter: Parameter already floating. No parameter changes made."
168  << std::endl;
169  }
170  } else {
171  std::cerr << "WARNING in LauPolarFormFactorSymNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
172  << std::endl;
173  }
174 }
175 
177 {
178  if ( name == "lambda" ) {
179  return lambda_;
180  } else {
181  std::cerr << "WARNING in LauPolarFormFactorSymNR::getResonanceParameter: Parameter name not reconised."
182  << std::endl;
183  return 0;
184  }
185 }
186 
187 void LauPolarFormFactorSymNR::setLambda( const Double_t lambda )
188 {
189  lambda_->value( lambda );
190  lambda_->genValue( lambda );
191  lambda_->initValue( lambda );
192 }
File containing declaration of LauResonanceInfo class.
const TString & getSanitisedName() const
Get the name of the resonance.
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.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
virtual ~LauPolarFormFactorSymNR()
Destructor.
File containing declaration of LauParameter class.
Double_t getLambda() const
Get the parameter lambda, the NR shape parameter.
File containing declaration of LauPolarFormFactorSymNR class.
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
Bool_t fixLambda() const
See if the lambda parameter is fixed or floating.
Int_t getSpin() const
Get the spin of the resonance.
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit.
Class for defining a complex number.
Definition: LauComplex.hh:61
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
const LauDaughters * getDaughters() const
Access the daughters object.
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 LauDaughters class.
void addExtraParameter(LauParameter *param, const Bool_t independentPar=kFALSE)
Add an extra parameter of the resonance.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
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.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
LauAbsResonance::LauResonanceModel model_
The model to use.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Get the complex dynamical amplitude.
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
Class for calculating 3-body kinematic quantities.
Double_t genValue() const
The value generated for the parameter.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getm23Sq() const
Get the m23 invariant mass square.
Double_t getm13Sq() const
Get the m13 invariant mass square.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:84
Double_t initValue() const
The initial value of the parameter.
LauResonanceModel
Define the allowed resonance types.