laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
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 "LauBelleSymNR.hh"
30 
31 #include "LauDaughters.hh"
32 #include "LauKinematics.hh"
33 #include "LauResonanceInfo.hh"
34 
35 #include "TMath.h"
36 
37 #include <iostream>
38 
41  const Int_t resPairAmpInt,
42  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."
66  << std::endl;
67  }
68 
69  Int_t resPairAmpInt = this->getPairInt();
70  if ( resPairAmpInt == 3 ) {
71  std::cerr << "WARNING in LauBelleSymNR::initialise : This lineshape is intended to be on the symmetrised axes of the DP."
72  << std::endl;
73  }
74 
75  if ( ( model_ != LauAbsResonance::BelleSymNR ) &&
78  std::cerr << "WARNING in LauBelleSymNR::initialise : Unknown model requested, defaulting to exponential."
79  << std::endl;
81  }
82 
83  if ( ( model_ != LauAbsResonance::BelleSymNRNoInter ) && ( this->getSpin() != 0 ) ) {
84  std::cerr << "WARNING in LauBelleSymNR::initialise : Non-zero spin will be ignored for this model - perhaps you should use LauAbsResonance::BelleSymNRNoInter instead"
85  << std::endl;
86  }
87 
88  // 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
89 }
90 
92 {
93  // This function returns the complex dynamical amplitude for a Belle Non-Resonant distribution
94 
95  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
96  // Have s<->t symmetry already done in Dynamics flip function.
97  // For Kpipi or similar plots, one can use the separate exponentials
98  // and consider them as two separate components with their own mag and phase.
99  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
100  // individual amplitudes (with the same value of alpha).
101 
102  // Calculate Mandelstam variables.
103  // s = m_13^2, t = m_23^2
104  const Double_t s = kinematics->getm13Sq();
105  const Double_t t = kinematics->getm23Sq();
106 
107  Double_t magnitude( 1.0 );
108 
109  const Double_t alpha = this->getAlpha();
110 
112 
113  magnitude = TMath::Exp( -alpha * s ) + TMath::Exp( -alpha * t );
114 
115  } else if ( model_ == LauAbsResonance::BelleSymNRNoInter ) {
116 
117  magnitude = ( s <= t ) ? TMath::Exp( -alpha * s ) : TMath::Exp( -alpha * t );
118 
119  const Double_t cosHel = ( s <= t ) ? kinematics->getc13() : kinematics->getc23();
120  magnitude *= this->calcLegendrePoly( cosHel );
121 
122  } else if ( model_ == LauAbsResonance::TaylorNR ) {
123 
124  const Double_t mParentSq = kinematics->getmParentSq();
125  magnitude = alpha * ( s + t ) / mParentSq + 1.0;
126  }
127 
128  LauComplex resAmplitude( magnitude, 0.0 );
129 
130  return resAmplitude;
131 }
132 
133 LauComplex LauBelleSymNR::resAmp( Double_t mass, Double_t spinTerm )
134 {
135  std::cerr << "ERROR in LauBelleSymNR : This method should never be called." << std::endl;
136  std::cerr << " : Returning zero amplitude for mass = " << mass
137  << " and spinTerm = " << spinTerm << "." << std::endl;
138  return LauComplex( 0.0, 0.0 );
139 }
140 
141 const std::vector<LauParameter*>& LauBelleSymNR::getFloatingParameters()
142 {
143  this->clearFloatingParameters();
144 
145  if ( ! this->fixAlpha() ) {
146  this->addFloatingParameter( alpha_ );
147  }
148 
149  return this->getParameters();
150 }
151 
152 void LauBelleSymNR::setResonanceParameter( const TString& name, const Double_t value )
153 {
154  // Set various parameters for the lineshape
155  if ( name == "alpha" ) {
156  this->setAlpha( value );
157  std::cout << "INFO in LauBelleSymNR::setResonanceParameter : Setting parameter alpha = "
158  << this->getAlpha() << std::endl;
159  } else {
160  std::cerr << "WARNING in LauBelleSymNR::setResonanceParameter : Parameter name not reconised. No parameter changes made."
161  << std::endl;
162  }
163 }
164 
166 {
167  if ( name == "alpha" ) {
168  if ( alpha_->fixed() ) {
169  alpha_->fixed( kFALSE );
170  this->addFloatingParameter( alpha_ );
171  } else {
172  std::cerr << "WARNING in LauBelleSymNR::floatResonanceParameter: Parameter already floating. No parameter changes made."
173  << std::endl;
174  }
175  } else {
176  std::cerr << "WARNING in LauBelleSymNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
177  << std::endl;
178  }
179 }
180 
182 {
183  if ( name == "alpha" ) {
184  return alpha_;
185  } else {
186  std::cerr << "WARNING in LauBelleSymNR::getResonanceParameter: Parameter name not reconised."
187  << std::endl;
188  return 0;
189  }
190 }
191 
192 void LauBelleSymNR::setAlpha( const Double_t alpha )
193 {
194  alpha_->value( alpha );
195  alpha_->genValue( alpha );
196  alpha_->initValue( alpha );
197 }
virtual void initialise()
Initialise.
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
Double_t getc13() const
Get the cosine of the helicity angle theta13.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
Double_t value() const
The value of the parameter.
LauAbsResonance::LauResonanceModel model_
The model to use.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
LauParameter * alpha_
The range parameter.
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
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.
Double_t getmParentSq() const
Get parent mass squared.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
const LauDaughters * getDaughters() const
Access the daughters object.
LauBelleSymNR(LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
virtual ~LauBelleSymNR()
Destructor.
File containing declaration of LauDaughters class.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
This is not called, amplitude is used directly instead.
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
File containing declaration of LauBelleSymNR class.
Class for defining the properties of a resonant particle.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
const TString & name() const
The parameter name.
Double_t getAlpha() const
Get the effective range parameter.
Bool_t fixAlpha() const
See if the alpha parameter is fixed or floating.
Double_t calcLegendrePoly() const
Calculate the Legendre polynomial for the spin factor.
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Class for calculating 3-body kinematic quantities.
Double_t genValue() const
The value generated for the parameter.
void setAlpha(const Double_t alpha)
Set the parameter alpha, the effective range.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getm23Sq() const
Get the m23 invariant mass square.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Double_t getm13Sq() const
Get the m13 invariant mass square.
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.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Get the complex dynamical amplitude.
File containing declaration of LauKinematics class.