laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauPolarFormFactorNR.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 
29 #include "LauPolarFormFactorNR.hh"
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  Int_t resPairAmpInt = this->getPairInt();
67  if ( daughters->gotSymmetricalDP() && resPairAmpInt != 3 ) {
68  std::cerr << "WARNING in LauPolarFormFactorNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate."
69  << std::endl;
70  }
71 
73  std::cerr << "WARNING in LauPolarFormFactorNR::initialise : Unknown model requested, defaulting to Polar Form Factor."
74  << std::endl;
76  }
77 }
78 
79 LauComplex LauPolarFormFactorNR::resAmp( Double_t mass, Double_t )
80 {
81  Double_t magnitude( 1.0 );
82 
83  Double_t lambda = this->getLambda();
84 
85  magnitude = 1.0 / ( 1.0 + mass * mass / ( lambda * lambda ) );
86 
87  LauComplex resAmplitude( magnitude, 0.0 );
88 
89  return resAmplitude;
90 }
91 
92 const std::vector<LauParameter*>& LauPolarFormFactorNR::getFloatingParameters()
93 {
95 
96  if ( ! this->fixLambda() ) {
98  }
99 
100  return this->getParameters();
101 }
102 
103 void LauPolarFormFactorNR::setResonanceParameter( const TString& name, const Double_t value )
104 {
105  // Set various parameters for the lineshape
106  if ( name == "lambda" ) {
107  this->setLambda( value );
108  std::cout << "INFO in LauPolarFormFactorNR::setResonanceParameter : Setting parameter lambda = "
109  << this->getLambda() << std::endl;
110  } else {
111  std::cerr << "WARNING in LauPolarFormFactorNR::setResonanceParameter: Parameter name not reconised. No parameter changes made."
112  << std::endl;
113  }
114 }
115 
117 {
118  if ( name == "lambda" ) {
119  if ( lambda_->fixed() ) {
120  lambda_->fixed( kFALSE );
121  this->addFloatingParameter( lambda_ );
122  } else {
123  std::cerr << "WARNING in LauPolarFormFactorNR::floatResonanceParameter: Parameter already floating. No parameter changes made."
124  << std::endl;
125  }
126  } else {
127  std::cerr << "WARNING in LauPolarFormFactorNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
128  << std::endl;
129  }
130 }
131 
133 {
134  if ( name == "lambda" ) {
135  return lambda_;
136  } else {
137  std::cerr << "WARNING in LauPolarFormFactorNR::getResonanceParameter: Parameter name not reconised."
138  << std::endl;
139  return 0;
140  }
141 }
142 
143 void LauPolarFormFactorNR::setLambda( const Double_t lambda )
144 {
145  lambda_->value( lambda );
146  lambda_->genValue( lambda );
147  lambda_->initValue( lambda );
148 }
File containing declaration of LauResonanceInfo class.
const TString & getSanitisedName() const
Get the name of the resonance.
Double_t getLambda() const
Get the parameter lambda, the NR shape parameter.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
Double_t value() const
The value of the parameter.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
virtual ~LauPolarFormFactorNR()
Destructor.
File containing declaration of LauParameter class.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
LauParameter * lambda_
The NR shape parameter.
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.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
const LauDaughters * getDaughters() const
Access the daughters object.
LauPolarFormFactorNR(LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
Bool_t fixLambda() const
See if the lambda parameter is fixed or floating.
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.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
void setLambda(const Double_t lambda)
Set the parameter lambda, the NR shape parameter.
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
Class for defining the properties of a resonant particle.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
virtual void initialise()
Initialise the model.
File containing declaration of LauPolarFormFactorNR class.
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Double_t genValue() const
The value generated for the parameter.
LauAbsResonance::LauResonanceModel model_
The model to use.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
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.