laura is hosted by Hepforge, IPPP Durham
Laura++  v3r5
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 <iostream>
30 
31 #include "TMath.h"
32 
33 #include "LauConstants.hh"
34 #include "LauPolarFormFactorNR.hh"
35 #include "LauDaughters.hh"
36 #include "LauParameter.hh"
37 #include "LauResonanceInfo.hh"
38 
40 
41 
42 LauPolarFormFactorNR::LauPolarFormFactorNR(LauResonanceInfo* resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters* daughters) :
43  LauAbsResonance(resInfo, resPairAmpInt, daughters),
44  lambda_(0),
45  model_(resType)
46 
47 {
48  TString parName = this->getSanitisedName();
49  parName += "_lambda";
50  lambda_ = resInfo->getExtraParameter( parName );
51  if ( lambda_ == 0 ) {
52  lambda_ = new LauParameter( parName, 1.0, 0.0, 10.0, kTRUE );
53  lambda_->secondStage(kTRUE);
54  resInfo->addExtraParameter( lambda_ );
55  }
56 }
57 
59 {
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." << std::endl;
69  }
70 
72  std::cerr << "WARNING in LauPolarFormFactorNR::initialise : Unknown model requested, defaulting to Polar Form Factor." << std::endl;
74  }
75 
76 }
77 
78 LauComplex LauPolarFormFactorNR::resAmp(Double_t mass, Double_t)
79 {
80  Double_t magnitude(1.0);
81 
82  Double_t lambda = this->getLambda();
83 
84  magnitude = 1.0/(1.0 + mass*mass /(lambda*lambda));
85 
86  LauComplex resAmplitude(magnitude, 0.0);
87 
88  return resAmplitude;
89 }
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 = " << this->getLambda() << std::endl;
109  }
110  else {
111  std::cerr << "WARNING in LauPolarFormFactorNR::setResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
112  }
113 }
114 
116 {
117  if (name == "lambda") {
118  if ( lambda_->fixed() ) {
119  lambda_->fixed( kFALSE );
120  this->addFloatingParameter( lambda_ );
121  } else {
122  std::cerr << "WARNING in LauPolarFormFactorNR::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl;
123  }
124  }
125  else {
126  std::cerr << "WARNING in LauPolarFormFactorNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
127  }
128 }
129 
131 {
132  if (name == "lambda") {
133  return lambda_;
134  }
135  else {
136  std::cerr << "WARNING in LauPolarFormFactorNR::getResonanceParameter: Parameter name not reconised." << std::endl;
137  return 0;
138  }
139 }
140 
141 void LauPolarFormFactorNR::setLambda(const Double_t lambda)
142 {
143  lambda_->value( lambda );
144  lambda_->genValue( lambda );
145  lambda_->initValue( lambda );
146 }
147 
148 
Bool_t fixed() const
Check whether the parameter is fixed or floated.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
File containing declaration of LauResonanceInfo class.
ClassImp(LauAbsCoeffSet)
LauParameter()
Default constructor.
Definition: LauParameter.cc:44
Class for defining the properties of a resonant particle.
const TString & name() const
The parameter name.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
File containing declaration of LauPolarFormFactorNR class.
File containing declaration of LauDaughters class.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
LauAbsResonance::LauResonanceModel model_
The model to use.
const LauDaughters * getDaughters() const
Access the daughters object.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
virtual ~LauPolarFormFactorNR()
Destructor.
Double_t getLambda() const
Get the parameter lambda, the NR shape parameter.
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:80
Class for defining a nonresonant form factor model.
File containing declaration of LauParameter class.
void setLambda(const Double_t lambda)
Set the parameter lambda, the NR shape parameter.
virtual void initialise()
Initialise the model.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
LauResonanceModel
Define the allowed resonance types.
Bool_t fixLambda() const
See if the lambda parameter is fixed or floating.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
Double_t initValue() const
The initial value of the parameter.
File containing LauConstants namespace.
Class for defining a complex number.
Definition: LauComplex.hh:61
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
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.
Double_t genValue() const
The value generated for the parameter.
void clearFloatingParameters()
Clear list of floating parameters.