laura is hosted by Hepforge, IPPP Durham
Laura++  v3r5
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauBelleNR.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2004 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 "LauBelleNR.hh"
34 #include "LauDaughters.hh"
35 #include "LauParameter.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  forceLegendre_(kTRUE)
47 {
48  TString parName = this->getSanitisedName();
49  parName += "_alpha";
50  alpha_ = resInfo->getExtraParameter( parName );
51  if ( alpha_ == 0 ) {
52  alpha_ = new LauParameter( parName, 0.0, 0.0, 10.0, kTRUE );
53  alpha_->secondStage(kTRUE);
54  resInfo->addExtraParameter( alpha_ );
55  }
56 }
57 
59 {
60 }
61 
63 {
64  const LauDaughters* daughters = this->getDaughters();
65  Int_t resPairAmpInt = this->getPairInt();
66  if ( daughters->gotSymmetricalDP() && resPairAmpInt != 3 ) {
67  std::cerr << "WARNING in LauBelleNR::initialise : Dalitz plot is symmetric - this lineshape is not appropriate." << std::endl;
68  }
69 
71  std::cerr << "WARNING in LauBelleNR::initialise : Unknown model requested, defaulting to exponential." << std::endl;
73  }
74 
75  // Make the spin term purely the Legendre polynomial of the cos(helicity angle)
76  if ( forceLegendre_ ) {
78  }
79 }
80 
81 LauComplex LauBelleNR::resAmp(Double_t mass, Double_t spinTerm)
82 {
83  Double_t magnitude(1.0);
84 
85  Double_t alpha = this->getAlpha();
86 
88  magnitude = spinTerm * TMath::Exp(-alpha*mass*mass);
89  } else if ( model_ == LauAbsResonance::PowerLawNR ) {
90  magnitude = spinTerm * TMath::Power(mass*mass, -alpha);
91  }
92 
93  LauComplex resAmplitude(magnitude, 0.0);
94 
95  return resAmplitude;
96 }
97 
98 const std::vector<LauParameter*>& LauBelleNR::getFloatingParameters()
99 {
100  this->clearFloatingParameters();
101 
102  if ( ! this->fixAlpha() ) {
103  this->addFloatingParameter( alpha_ );
104  }
105 
106  return this->getParameters();
107 }
108 
109 void LauBelleNR::setResonanceParameter(const TString& name, const Double_t value)
110 {
111  // Set various parameters for the lineshape
112  if (name == "alpha") {
113  this->setAlpha(value);
114  std::cout << "INFO in LauBelleNR::setResonanceParameter : Setting parameter alpha = " << this->getAlpha() << std::endl;
115  } else {
116  std::cerr << "WARNING in LauBelleNR::setResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
117  }
118 }
119 
121 {
122  if (name == "alpha") {
123  if ( alpha_->fixed() ) {
124  alpha_->fixed( kFALSE );
125  this->addFloatingParameter( alpha_ );
126  } else {
127  std::cerr << "WARNING in LauBelleNR::floatResonanceParameter: Parameter already floating. No parameter changes made." << std::endl;
128  }
129  } else {
130  std::cerr << "WARNING in LauBelleNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made." << std::endl;
131  }
132 }
133 
135 {
136  if (name == "alpha") {
137  return alpha_;
138  } else {
139  std::cerr << "WARNING in LauBelleNR::getResonanceParameter: Parameter name not reconised." << std::endl;
140  return 0;
141  }
142 }
143 
144 void LauBelleNR::setAlpha(const Double_t alpha)
145 {
146  alpha_->value( alpha );
147  alpha_->genValue( alpha );
148  alpha_->initValue( alpha );
149 }
150 
LauParameter * alpha_
The range parameter.
Definition: LauBelleNR.hh:144
Bool_t fixed() const
Check whether the parameter is fixed or floated.
virtual ~LauBelleNR()
Destructor.
Definition: LauBelleNR.cc:58
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
File containing declaration of LauBelleNR class.
File containing declaration of LauDaughters class.
Bool_t forceLegendre_
Force use of Legendre spin factors.
Definition: LauBelleNR.hh:150
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
Definition: LauBelleNR.cc:120
Double_t getAlpha() const
Get the effective range parameter.
Definition: LauBelleNR.hh:121
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
Definition: LauBelleNR.cc:81
const LauDaughters * getDaughters() const
Access the daughters object.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Definition: LauBelleNR.cc:98
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:80
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Definition: LauBelleNR.cc:109
virtual void initialise()
Initialise the model.
Definition: LauBelleNR.cc:62
Class for defining the Belle nonresonant model.
Definition: LauBelleNR.hh:48
File containing declaration of LauParameter class.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
Bool_t fixAlpha() const
See if the alpha parameter is fixed or floating.
Definition: LauBelleNR.hh:127
LauResonanceModel
Define the allowed resonance types.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Definition: LauBelleNR.cc:134
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
void setAlpha(const Double_t alpha)
Set the parameter alpha, the effective range.
Definition: LauBelleNR.cc:144
Double_t initValue() const
The initial value of the parameter.
Class for defining a complex number.
Definition: LauComplex.hh:61
LauAbsResonance::LauResonanceModel model_
The model to use.
Definition: LauBelleNR.hh:147
void setSpinType(const LauSpinType spinType)
Set the spin formalism to be used.
Double_t value() const
The value of the parameter.
Double_t genValue() const
The value generated for the parameter.
void clearFloatingParameters()
Clear list of floating parameters.