laura is hosted by Hepforge, IPPP Durham
Laura++  v3r2
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauBelleSymNR.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2014.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 
17 #include "TMath.h"
18 
19 #include "LauBelleSymNR.hh"
20 #include "LauDaughters.hh"
21 #include "LauKinematics.hh"
22 #include "LauResonanceInfo.hh"
23 
25 
26 
28  const Int_t resPairAmpInt, const LauDaughters* daughters) :
29  LauAbsResonance(resInfo, resPairAmpInt, daughters),
30  alpha_(0),
31  model_(resType)
32 {
33  TString parName = this->getSanitisedName();
34  parName += "_alpha";
35  alpha_ = resInfo->getExtraParameter( parName );
36  if ( alpha_ == 0 ) {
37  alpha_ = new LauParameter( parName, 0.0, 0.0, 10.0, kTRUE );
38  alpha_->secondStage(kTRUE);
39  resInfo->addExtraParameter( alpha_ );
40  }
41 }
42 
44 {
45 }
46 
48 {
49  const LauDaughters* daughters = this->getDaughters();
50  if ( ! daughters->gotSymmetricalDP() ) {
51  std::cerr << "WARNING in LauBelleSymNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate." << std::endl;
52  }
53 
54  Int_t resPairAmpInt = this->getPairInt();
55  if ( resPairAmpInt == 3 ) {
56  std::cerr << "WARNING in LauBelleSymNR::initialise : This lineshape is intended to be on the symmetrised axes of the DP." << std::endl;
57  }
58 
60  std::cerr << "WARNING in LauBelleSymNR::initialise : Unknown model requested, defaulting to exponential." << std::endl;
62  }
63 
64  if ( (model_ != LauAbsResonance::BelleSymNRNoInter) && (this->getSpin() != 0) ) {
65  std::cerr << "WARNING in LauBelleSymNR::initialise : Non-zero spin will be ignored for this model - perhaps you should use LauAbsResonance::BelleSymNRNoInter instead" << std::endl;
66  }
67 
68  // 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
69 }
70 
72 {
73  // This function returns the complex dynamical amplitude for a Belle Non-Resonant distribution
74 
75  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
76  // Have s<->t symmetry already done in Dynamics flip function.
77  // For Kpipi or similar plots, one can use the separate exponentials
78  // and consider them as two separate components with their own mag and phase.
79  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
80  // individual amplitudes (with the same value of alpha).
81 
82  // Calculate Mandelstam variables.
83  // s = m_13^2, t = m_23^2
84  const Double_t s = kinematics->getm13Sq();
85  const Double_t t = kinematics->getm23Sq();
86 
87  Double_t magnitude(1.0);
88 
89  const Double_t alpha = this->getAlpha();
90 
92 
93  magnitude = TMath::Exp(-alpha*s) + TMath::Exp(-alpha*t);
94 
96 
97  magnitude = (s <= t) ? TMath::Exp(-alpha*s) : TMath::Exp(-alpha*t);
98 
99  const Double_t cosHel = (s <= t) ? kinematics->getc13() : kinematics->getc23();
100  magnitude *= this->calcLegendrePoly( cosHel );
101 
102  } else if ( model_ == LauAbsResonance::TaylorNR ) {
103 
104  const Double_t mParentSq = kinematics->getmParentSq();
105  magnitude = alpha*(s + t)/mParentSq + 1.0;
106 
107  }
108 
109  LauComplex resAmplitude(magnitude, 0.0);
110 
111  return resAmplitude;
112 }
113 
114 LauComplex LauBelleSymNR::resAmp(Double_t mass, Double_t spinTerm)
115 {
116  std::cerr << "ERROR in LauBelleSymNR : This method should never be called." << std::endl;
117  std::cerr << " : Returning zero amplitude for mass = " << mass << " and spinTerm = " << spinTerm << "." << std::endl;
118  return LauComplex(0.0, 0.0);
119 }
120 
121 const std::vector<LauParameter*>& LauBelleSymNR::getFloatingParameters()
122 {
123  this->clearFloatingParameters();
124 
125  if ( ! this->fixAlpha() ) {
126  this->addFloatingParameter( alpha_ );
127  }
128 
129  return this->getParameters();
130 }
131 
132 void LauBelleSymNR::setResonanceParameter(const TString& name, const Double_t value)
133 {
134  // Set various parameters for the lineshape
135  if (name == "alpha") {
136  this->setAlpha(value);
137  std::cout << "INFO in LauBelleSymNR::setResonanceParameter : Setting parameter alpha = " << this->getAlpha() << std::endl;
138  } else {
139  std::cerr << "WARNING in LauBelleSymNR::setResonanceParameter : Parameter name not reconised. No parameter changes made." << std::endl;
140  }
141 }
142 
144 {
145  if (name == "alpha") {
146  if ( alpha_->fixed() ) {
147  alpha_->fixed( kFALSE );
148  this->addFloatingParameter( alpha_ );
149  } else {
150  std::cerr << "WARNING in LauBelleSymNR::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl;
151  }
152  } else {
153  std::cerr << "WARNING in LauBelleSymNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
154  }
155 }
156 
158 {
159  if (name == "alpha") {
160  return alpha_;
161  } else {
162  std::cerr << "WARNING in LauBelleSymNR::getResonanceParameter: Parameter name not reconised." << std::endl;
163  return 0;
164  }
165 }
166 
167 void LauBelleSymNR::setAlpha(const Double_t alpha)
168 {
169  alpha_->value( alpha );
170  alpha_->genValue( alpha );
171  alpha_->initValue( alpha );
172 }
173 
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:30
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:33
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:66
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:35
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:47
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.