laura is hosted by Hepforge, IPPP Durham
Laura++  v2r2p1
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 - 2013.
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 
24 
25 
26 LauBelleSymNR::LauBelleSymNR(const TString& resName, const LauAbsResonance::LauResonanceModel resType,
27  const Double_t resMass, const Double_t resWidth, const Int_t resSpin,
28  const Int_t resCharge, const Int_t resPairAmpInt, const LauDaughters* daughters) :
29  LauAbsResonance(resName, resMass, resWidth, resSpin, resCharge, resPairAmpInt, daughters),
30  alpha_(0.0),
31  model_(resType)
32 {
33 }
34 
36 {
37 }
38 
40 {
41  const LauDaughters* daughters = this->getDaughters();
42  Int_t resPairAmpInt = this->getPairInt();
43  if ( ! daughters->gotSymmetricalDP() ) {
44  std::cerr << "WARNING in LauBelleSymNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate." << std::endl;
45  }
46  if ( resPairAmpInt == 3 ) {
47  std::cerr << "WARNING in LauBelleSymNR::initialise : This lineshape is intended to be on the symmetrised axes of the DP." << std::endl;
48  }
49 
51  std::cerr << "WARNING in LauBelleSymNR::initialise : Unknown model requested, defaulting to exponential." << std::endl;
53  }
54 }
55 
57 {
58  // This function returns the complex dynamical amplitude for a Belle Non-Resonant distribution
59 
60  // Calculate for symmetric DPs, e.g. 3pi or 3K, by using shapeNo = 1 or 2
61  // Have s<->t symmetry already done in Dynamics flip function.
62  // For Kpipi or similar plots, one can use the separate exponentials
63  // and consider them as two separate components with their own mag and phase.
64  // For this shapeNo = 3 and shapeNo = 4 need to be used to create the two
65  // individual amplitudes (with the same value of alpha).
66 
67  // Calculate Mandelstam variables.
68  // s = m_13^2, t = m_23^2, u = m_12^2.
69  Double_t s = kinematics->getm13Sq();
70  Double_t t = kinematics->getm23Sq();
71  //Double_t u = kinematics->getm12Sq();
72 
73  Double_t magnitude(1.0);
74 
76  magnitude = TMath::Exp(-alpha_*s) + TMath::Exp(-alpha_*t);
77  } else if ( model_ == LauAbsResonance::TaylorNR ) {
78  Double_t mParentSq = kinematics->getmParentSq();
79  magnitude = alpha_*(s + t)/mParentSq + 1.0;
80  }
81 
82  LauComplex resAmplitude(magnitude, 0.0);
83 
84  return resAmplitude;
85 }
86 
87 LauComplex LauBelleSymNR::resAmp(Double_t mass, Double_t spinTerm)
88 {
89  std::cerr << "ERROR in LauBelleSymNR : This method should never be called." << std::endl;
90  std::cerr << " : Returning zero amplitude for mass = " << mass << " and spinTerm = " << spinTerm << "." << std::endl;
91  return LauComplex(0.0, 0.0);
92 }
93 
94 void LauBelleSymNR::setResonanceParameter(const TString& name, const Double_t value)
95 {
96  // Set various parameters for the lineshape
97  if (name == "alpha") {
98  this->setAlpha(value);
99  std::cout << "INFO in LauBelleSymNR::setResonanceParameter : Setting parameter alpha = " << this->getAlpha() << std::endl;
100  } else {
101  std::cerr << "WARNING in LauBelleSymNR::setResonanceParameter : Parameter name not reconised. No parameter changes made." << std::endl;
102  }
103 }
104 
virtual void setAlpha(Double_t alpha)
Set the parameter alpha, the effective range.
ClassImp(LauAbsCoeffSet)
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.
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.
virtual ~LauBelleSymNR()
Destructor.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric.
Definition: LauDaughters.hh:66
virtual Double_t getAlpha()
Get the effective range parameter.
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.
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 getmParentSq() const
Get parent mass squared.
Double_t alpha_
Define the range parameter.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.