laura is hosted by Hepforge, IPPP Durham
Laura++  v3r0
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  delete alpha_;
46 }
47 
49 {
50  const LauDaughters* daughters = this->getDaughters();
51  Int_t resPairAmpInt = this->getPairInt();
52  if ( ! daughters->gotSymmetricalDP() ) {
53  std::cerr << "WARNING in LauBelleSymNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate." << std::endl;
54  }
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 
66 {
67  // This function returns the complex dynamical amplitude for a Belle Non-Resonant distribution
68 
69  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
70  // Have s<->t symmetry already done in Dynamics flip function.
71  // For Kpipi or similar plots, one can use the separate exponentials
72  // and consider them as two separate components with their own mag and phase.
73  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
74  // individual amplitudes (with the same value of alpha).
75 
76  // Calculate Mandelstam variables.
77  // s = m_13^2, t = m_23^2, u = m_12^2.
78  Double_t s = kinematics->getm13Sq();
79  Double_t t = kinematics->getm23Sq();
80  //Double_t u = kinematics->getm12Sq();
81 
82  Double_t magnitude(1.0);
83 
84  Double_t alpha = this->getAlpha();
85 
87  magnitude = TMath::Exp(-alpha*s) + TMath::Exp(-alpha*t);
88  } else if ( model_ == LauAbsResonance::TaylorNR ) {
89  Double_t mParentSq = kinematics->getmParentSq();
90  magnitude = alpha*(s + t)/mParentSq + 1.0;
91  }
92 
93  LauComplex resAmplitude(magnitude, 0.0);
94 
95  return resAmplitude;
96 }
97 
98 LauComplex LauBelleSymNR::resAmp(Double_t mass, Double_t spinTerm)
99 {
100  std::cerr << "ERROR in LauBelleSymNR : This method should never be called." << std::endl;
101  std::cerr << " : Returning zero amplitude for mass = " << mass << " and spinTerm = " << spinTerm << "." << std::endl;
102  return LauComplex(0.0, 0.0);
103 }
104 
105 const std::vector<LauParameter*>& LauBelleSymNR::getFloatingParameters()
106 {
107  this->clearFloatingParameters();
108 
109  if ( ! this->fixAlpha() ) {
110  this->addFloatingParameter( alpha_ );
111  }
112 
113  return this->getParameters();
114 }
115 
116 void LauBelleSymNR::setResonanceParameter(const TString& name, const Double_t value)
117 {
118  // Set various parameters for the lineshape
119  if (name == "alpha") {
120  this->setAlpha(value);
121  std::cout << "INFO in LauBelleSymNR::setResonanceParameter : Setting parameter alpha = " << this->getAlpha() << std::endl;
122  } else {
123  std::cerr << "WARNING in LauBelleSymNR::setResonanceParameter : Parameter name not reconised. No parameter changes made." << std::endl;
124  }
125 }
126 
128 {
129  if (name == "alpha") {
130  if ( alpha_->fixed() ) {
131  alpha_->fixed( kFALSE );
132  this->addFloatingParameter( alpha_ );
133  } else {
134  std::cerr << "WARNING in LauBelleSymNR::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl;
135  }
136  } else {
137  std::cerr << "WARNING in LauBelleSymNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
138  }
139 }
140 
142 {
143  if (name == "alpha") {
144  return alpha_;
145  } else {
146  std::cerr << "WARNING in LauBelleSymNR::getResonanceParameter: Parameter name not reconised." << std::endl;
147  return 0;
148  }
149 }
150 
151 void LauBelleSymNR::setAlpha(const Double_t alpha)
152 {
153  alpha_->value( alpha );
154  alpha_->genValue( alpha );
155  alpha_->initValue( alpha );
156 }
157 
Bool_t fixed() const
Check whether the parameter is fixed or floated.
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
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:34
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.
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.