laura is hosted by Hepforge, IPPP Durham
Laura++  v3r5
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauBelleSymNR.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2013 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 "LauBelleSymNR.hh"
34 #include "LauDaughters.hh"
35 #include "LauKinematics.hh"
36 #include "LauResonanceInfo.hh"
37 
39 
40 
42  const Int_t resPairAmpInt, const LauDaughters* daughters) :
43  LauAbsResonance(resInfo, resPairAmpInt, daughters),
44  alpha_(0),
45  model_(resType)
46 {
47  TString parName = this->getSanitisedName();
48  parName += "_alpha";
49  alpha_ = resInfo->getExtraParameter( parName );
50  if ( alpha_ == 0 ) {
51  alpha_ = new LauParameter( parName, 0.0, 0.0, 10.0, kTRUE );
52  alpha_->secondStage(kTRUE);
53  resInfo->addExtraParameter( alpha_ );
54  }
55 }
56 
58 {
59 }
60 
62 {
63  const LauDaughters* daughters = this->getDaughters();
64  if ( ! daughters->gotSymmetricalDP() ) {
65  std::cerr << "WARNING in LauBelleSymNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate." << std::endl;
66  }
67 
68  Int_t resPairAmpInt = this->getPairInt();
69  if ( resPairAmpInt == 3 ) {
70  std::cerr << "WARNING in LauBelleSymNR::initialise : This lineshape is intended to be on the symmetrised axes of the DP." << std::endl;
71  }
72 
74  std::cerr << "WARNING in LauBelleSymNR::initialise : Unknown model requested, defaulting to exponential." << std::endl;
76  }
77 
78  if ( (model_ != LauAbsResonance::BelleSymNRNoInter) && (this->getSpin() != 0) ) {
79  std::cerr << "WARNING in LauBelleSymNR::initialise : Non-zero spin will be ignored for this model - perhaps you should use LauAbsResonance::BelleSymNRNoInter instead" << std::endl;
80  }
81 
82  // NB we do not need to call setSpinType(LauAbsResonance::Legendre) here (as is done in LauBelleNR) since override the amplitude method and explicitly use calcLegendrePoly
83 }
84 
86 {
87  // This function returns the complex dynamical amplitude for a Belle Non-Resonant distribution
88 
89  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
90  // Have s<->t symmetry already done in Dynamics flip function.
91  // For Kpipi or similar plots, one can use the separate exponentials
92  // and consider them as two separate components with their own mag and phase.
93  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
94  // individual amplitudes (with the same value of alpha).
95 
96  // Calculate Mandelstam variables.
97  // s = m_13^2, t = m_23^2
98  const Double_t s = kinematics->getm13Sq();
99  const Double_t t = kinematics->getm23Sq();
100 
101  Double_t magnitude(1.0);
102 
103  const Double_t alpha = this->getAlpha();
104 
106 
107  magnitude = TMath::Exp(-alpha*s) + TMath::Exp(-alpha*t);
108 
109  } else if ( model_ == LauAbsResonance::BelleSymNRNoInter ) {
110 
111  magnitude = (s <= t) ? TMath::Exp(-alpha*s) : TMath::Exp(-alpha*t);
112 
113  const Double_t cosHel = (s <= t) ? kinematics->getc13() : kinematics->getc23();
114  magnitude *= this->calcLegendrePoly( cosHel );
115 
116  } else if ( model_ == LauAbsResonance::TaylorNR ) {
117 
118  const Double_t mParentSq = kinematics->getmParentSq();
119  magnitude = alpha*(s + t)/mParentSq + 1.0;
120 
121  }
122 
123  LauComplex resAmplitude(magnitude, 0.0);
124 
125  return resAmplitude;
126 }
127 
128 LauComplex LauBelleSymNR::resAmp(Double_t mass, Double_t spinTerm)
129 {
130  std::cerr << "ERROR in LauBelleSymNR : This method should never be called." << std::endl;
131  std::cerr << " : Returning zero amplitude for mass = " << mass << " and spinTerm = " << spinTerm << "." << std::endl;
132  return LauComplex(0.0, 0.0);
133 }
134 
135 const std::vector<LauParameter*>& LauBelleSymNR::getFloatingParameters()
136 {
137  this->clearFloatingParameters();
138 
139  if ( ! this->fixAlpha() ) {
140  this->addFloatingParameter( alpha_ );
141  }
142 
143  return this->getParameters();
144 }
145 
146 void LauBelleSymNR::setResonanceParameter(const TString& name, const Double_t value)
147 {
148  // Set various parameters for the lineshape
149  if (name == "alpha") {
150  this->setAlpha(value);
151  std::cout << "INFO in LauBelleSymNR::setResonanceParameter : Setting parameter alpha = " << this->getAlpha() << std::endl;
152  } else {
153  std::cerr << "WARNING in LauBelleSymNR::setResonanceParameter : Parameter name not reconised. No parameter changes made." << std::endl;
154  }
155 }
156 
158 {
159  if (name == "alpha") {
160  if ( alpha_->fixed() ) {
161  alpha_->fixed( kFALSE );
162  this->addFloatingParameter( alpha_ );
163  } else {
164  std::cerr << "WARNING in LauBelleSymNR::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl;
165  }
166  } else {
167  std::cerr << "WARNING in LauBelleSymNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
168  }
169 }
170 
172 {
173  if (name == "alpha") {
174  return alpha_;
175  } else {
176  std::cerr << "WARNING in LauBelleSymNR::getResonanceParameter: Parameter name not reconised." << std::endl;
177  return 0;
178  }
179 }
180 
181 void LauBelleSymNR::setAlpha(const Double_t alpha)
182 {
183  alpha_->value( alpha );
184  alpha_->genValue( alpha );
185  alpha_->initValue( alpha );
186 }
187 
Double_t getc23() const
Get the cosine of the helicity angle theta23.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
Double_t calcLegendrePoly() const
Calculate the Legendre polynomial for the spin factor.
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
Double_t getc13() const
Get the cosine of the helicity angle theta13.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
This is not called, amplitude is used directly instead.
File containing declaration of LauDaughters class.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
File containing declaration of LauBelleSymNR class.
Double_t getm23Sq() const
Get the m23 invariant mass square.
const LauDaughters * getDaughters() const
Access the daughters object.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
File containing declaration of LauKinematics class.
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
virtual ~LauBelleSymNR()
Destructor.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:80
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
Double_t getAlpha() const
Get the effective range parameter.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
virtual LauComplex amplitude(const LauKinematics *kinematics)
Get the complex dynamical amplitude.
LauResonanceModel
Define the allowed resonance types.
Class for defining the symmetric Belle Non Resonant model.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
Double_t getm13Sq() const
Get the m13 invariant mass square.
Double_t initValue() const
The initial value of the parameter.
LauParameter * alpha_
The range parameter.
Bool_t fixAlpha() const
See if the alpha parameter is fixed or floating.
Class for defining a complex number.
Definition: LauComplex.hh:61
Class for calculating 3-body kinematic quantities.
Double_t value() const
The value of the parameter.
virtual void initialise()
Initialise.
LauAbsResonance::LauResonanceModel model_
The model to use.
Int_t getSpin() const
Get the spin of the resonance.
Double_t genValue() const
The value generated for the parameter.
Double_t getmParentSq() const
Get parent mass squared.
void setAlpha(const Double_t alpha)
Set the parameter alpha, the effective range.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
void clearFloatingParameters()
Clear list of floating parameters.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.