laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
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 "LauBelleNR.hh"
30 
31 #include "LauDaughters.hh"
32 #include "LauParameter.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  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, -2.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."
68  << std::endl;
69  }
70 
72  std::cerr << "WARNING in LauBelleNR::initialise : Unknown model requested, defaulting to exponential."
73  << std::endl;
75  }
76 
77  // Make the spin term purely the Legendre polynomial of the cos(helicity angle)
78  if ( forceLegendre_ ) {
80  }
81 }
82 
83 LauComplex LauBelleNR::resAmp( Double_t mass, Double_t spinTerm )
84 {
85  Double_t magnitude( 1.0 );
86 
87  Double_t alpha = this->getAlpha();
88 
90  magnitude = spinTerm * TMath::Exp( -alpha * mass * mass );
91  } else if ( model_ == LauAbsResonance::PowerLawNR ) {
92  magnitude = spinTerm * TMath::Power( mass * mass, -alpha );
93  }
94 
95  LauComplex resAmplitude( magnitude, 0.0 );
96 
97  return resAmplitude;
98 }
99 
100 const std::vector<LauParameter*>& LauBelleNR::getFloatingParameters()
101 {
102  this->clearFloatingParameters();
103 
104  if ( ! this->fixAlpha() ) {
105  this->addFloatingParameter( alpha_ );
106  }
107 
108  return this->getParameters();
109 }
110 
111 void LauBelleNR::setResonanceParameter( const TString& name, const Double_t value )
112 {
113  // Set various parameters for the lineshape
114  if ( name == "alpha" ) {
115  this->setAlpha( value );
116  std::cout << "INFO in LauBelleNR::setResonanceParameter : Setting parameter alpha = "
117  << this->getAlpha() << std::endl;
118  } else {
119  std::cerr << "WARNING in LauBelleNR::setResonanceParameter: Parameter name not reconised. No parameter changes made."
120  << std::endl;
121  }
122 }
123 
125 {
126  if ( name == "alpha" ) {
127  if ( alpha_->fixed() ) {
128  alpha_->fixed( kFALSE );
129  this->addFloatingParameter( alpha_ );
130  } else {
131  std::cerr << "WARNING in LauBelleNR::floatResonanceParameter: Parameter already floating. No parameter changes made."
132  << std::endl;
133  }
134  } else {
135  std::cerr << "WARNING in LauBelleNR::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
136  << std::endl;
137  }
138 }
139 
141 {
142  if ( name == "alpha" ) {
143  return alpha_;
144  } else {
145  std::cerr << "WARNING in LauBelleNR::getResonanceParameter: Parameter name not reconised."
146  << std::endl;
147  return 0;
148  }
149 }
150 
151 void LauBelleNR::setAlpha( const Double_t alpha )
152 {
153  alpha_->value( alpha );
154  alpha_->genValue( alpha );
155  alpha_->initValue( alpha );
156 }
virtual void initialise()
Initialise the model.
Definition: LauBelleNR.cc:62
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 value() const
The value of the parameter.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Definition: LauBelleNR.cc:140
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
File containing declaration of LauParameter class.
File containing declaration of LauBelleNR class.
Double_t getAlpha() const
Get the effective range parameter.
Definition: LauBelleNR.hh:125
void setAlpha(const Double_t alpha)
Set the parameter alpha, the effective range.
Definition: LauBelleNR.cc:151
LauAbsResonance::LauResonanceModel model_
The model to use.
Definition: LauBelleNR.hh:151
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter 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.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
const LauDaughters * getDaughters() const
Access the daughters object.
LauBelleNR(LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
Definition: LauBelleNR.cc:39
LauParameter * alpha_
The range parameter.
Definition: LauBelleNR.hh:148
File containing declaration of LauDaughters class.
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.
virtual ~LauBelleNR()
Destructor.
Definition: LauBelleNR.cc:58
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
Class for defining the properties of a resonant particle.
const TString & name() const
The parameter name.
Bool_t fixAlpha() const
See if the alpha parameter is fixed or floating.
Definition: LauBelleNR.hh:131
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
Definition: LauBelleNR.cc:83
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
Definition: LauBelleNR.cc:124
Bool_t forceLegendre_
Force use of Legendre spin factors.
Definition: LauBelleNR.hh:154
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Definition: LauBelleNR.cc:100
Double_t genValue() const
The value generated for the parameter.
void setSpinType(const LauSpinType spinType)
Set the spin formalism to be used.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Definition: LauBelleNR.cc:111
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.