laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauLASSRes.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 "LauLASSRes.hh"
30 
31 #include "LauConstants.hh"
32 #include "LauResonanceInfo.hh"
33 
34 #include <iostream>
35 
37  const Int_t resPairAmpInt,
38  const LauDaughters* daughters ) :
39  LauAbsResonance( resInfo, resPairAmpInt, daughters ),
40  q0_( 0.0 ),
41  mDaugSum_( 0.0 ),
42  mDaugSumSq_( 0.0 ),
43  mDaugDiff_( 0.0 ),
44  mDaugDiffSq_( 0.0 ),
45  resMassSq_( 0.0 ),
46  r_( 0 ),
47  a_( 0 ),
48  cutOff_( 0.0 )
49 {
50  // Default values for LASS parameters
51  cutOff_ = 1.8;
52  const Double_t rVal = 3.32;
53  const Double_t aVal = 2.07;
54 
55  const TString& parNameBase = this->getSanitisedName();
56 
57  TString rName( parNameBase );
58  rName += "_r";
59  r_ = resInfo->getExtraParameter( rName );
60  if ( r_ == 0 ) {
61  r_ = new LauParameter( rName, rVal, 0.0, 10.0, kTRUE );
62  r_->secondStage( kTRUE );
63  resInfo->addExtraParameter( r_ );
64  }
65 
66  TString aName( parNameBase );
67  aName += "_a";
68  a_ = resInfo->getExtraParameter( aName );
69  if ( a_ == 0 ) {
70  a_ = new LauParameter( aName, aVal, 0.0, 10.0, kTRUE );
71  a_->secondStage( kTRUE );
72  resInfo->addExtraParameter( a_ );
73  }
74 }
75 
77 {
78 }
79 
81 {
82  // Create the mass sums and differences
83  Double_t massDaug1 = this->getMassDaug1();
84  Double_t massDaug2 = this->getMassDaug2();
85 
86  mDaugSum_ = massDaug1 + massDaug2;
88 
89  mDaugDiff_ = massDaug1 - massDaug2;
91 
92  Int_t resSpin = this->getSpin();
93  if ( resSpin != 0 ) {
94  std::cerr << "WARNING in LauLASSRes::amplitude : Resonance spin is " << resSpin << "."
95  << std::endl;
96  std::cerr << " : LASS amplitude is only for scalers, resetting spin to 0."
97  << std::endl;
98  this->changeResonance( -1.0, -1.0, 0 );
99  }
100 
101  this->calcQ0();
102 }
103 
105 {
106  // Decay momentum of either daughter in the resonance rest frame
107  // when resonance mass = rest-mass value, m_0 (PDG value)
108 
109  resMass_ = this->getMass();
111 
112  q0_ = TMath::Sqrt( ( resMassSq_ - mDaugSumSq_ ) * ( resMassSq_ - mDaugDiffSq_ ) ) /
113  ( 2.0 * resMass_ );
114 }
115 
116 LauComplex LauLASSRes::resAmp( Double_t mass, Double_t /*spinTerm*/ )
117 {
118  LauComplex resAmplitude( 0.0, 0.0 );
119  LauComplex bkgAmplitude( 0.0, 0.0 );
120  LauComplex totAmplitude( 0.0, 0.0 );
121 
122  if ( mass < 1e-10 ) {
123  std::cerr << "WARNING in LauLASSRes::amplitude : mass < 1e-10." << std::endl;
124  return LauComplex( 0.0, 0.0 );
125  }
126 
127  //---------------------------
128  // First do the resonant part
129  //---------------------------
130 
131  // Calculate the width of the resonance (as a function of mass)
132  // q is the momentum of either daughter in the resonance rest-frame
133  const Double_t q = this->getQ();
134  const Double_t resMass = this->getMass();
135  const Double_t resWidth = this->getWidth();
136 
137  // If the mass is floating and their value have changed
138  // we need to recalculate everything that assumes this value
139  if ( ( ! this->fixMass() ) && resMass != resMass_ ) {
140  this->calcQ0();
141  }
142 
143  Double_t qRatio = q / q0_;
144  Double_t totWidth = resWidth * qRatio * ( resMass / mass );
145 
146  Double_t massSq = mass * mass;
147  Double_t massSqTerm = resMassSq_ - massSq;
148 
149  // Compute the complex amplitude
150  resAmplitude = LauComplex( massSqTerm, resMass * totWidth );
151 
152  // Scale by the denominator factor
153  resAmplitude.rescale( ( resMassSq_ * resWidth / q0_ ) /
154  ( massSqTerm * massSqTerm + resMassSq_ * totWidth * totWidth ) );
155 
156  // Calculate the phase shift term
157  const Double_t rVal = this->getEffectiveRange();
158  const Double_t aVal = this->getScatteringLength();
159  const Double_t tandeltaB = ( 2.0 * aVal * q ) / ( 2.0 + aVal * rVal * q * q );
160  const Double_t tanSq = tandeltaB * tandeltaB;
161  const Double_t cos2PhaseShift = ( 1.0 - tanSq ) / ( 1.0 + tanSq );
162  const Double_t sin2PhaseShift = 2.0 * tandeltaB / ( 1.0 + tanSq );
163  LauComplex phaseShift( cos2PhaseShift, sin2PhaseShift );
164 
165  // Multiply by the phase shift
166  resAmplitude = resAmplitude * phaseShift;
167 
168  //--------------------------------
169  // Now do the effective range part
170  //--------------------------------
171 
172  // Form the real and imaginary parts
173  const Double_t qcotdeltaB = 1.0 / aVal + ( rVal * q * q ) / 2.0;
174 
175  // Compute the complex amplitude
176  bkgAmplitude = LauComplex( qcotdeltaB, q );
177 
178  // Scale by the numerator and denominator factors
179  bkgAmplitude.rescale( mass / ( qcotdeltaB * qcotdeltaB + q * q ) );
180 
181  //------------------
182  // Add them together
183  //------------------
184 
185  if ( mass > cutOff_ ) {
186  totAmplitude = resAmplitude;
187  } else {
188  totAmplitude = bkgAmplitude + resAmplitude;
189  }
190 
191  return totAmplitude;
192 }
193 
194 const std::vector<LauParameter*>& LauLASSRes::getFloatingParameters()
195 {
196  this->clearFloatingParameters();
197 
198  if ( ! this->fixMass() ) {
199  this->addFloatingParameter( this->getMassPar() );
200  }
201 
202  if ( ! this->fixWidth() ) {
203  this->addFloatingParameter( this->getWidthPar() );
204  }
205 
206  if ( ! this->fixEffectiveRange() ) {
207  this->addFloatingParameter( r_ );
208  }
209 
210  if ( ! this->fixScatteringLength() ) {
211  this->addFloatingParameter( a_ );
212  }
213 
214  return this->getParameters();
215 }
216 
217 void LauLASSRes::setResonanceParameter( const TString& name, const Double_t value )
218 {
219  // Set various parameters for the LASS lineshape dynamics
220  if ( name == "a" ) {
221  this->setScatteringLength( value );
222  std::cout << "INFO in LauLASSRes::setResonanceParameter : Setting LASS Scattering Length = "
223  << this->getScatteringLength() << std::endl;
224  } else if ( name == "r" ) {
225  this->setEffectiveRange( value );
226  std::cout << "INFO in LauLASSRes::setResonanceParameter : Setting LASS Effective Range = "
227  << this->getEffectiveRange() << std::endl;
228  } else {
229  std::cerr << "WARNING in LauLASSRes::setResonanceParameter: Parameter name not reconised. No parameter changes made."
230  << std::endl;
231  }
232 }
233 
235 {
236  if ( name == "a" ) {
237  if ( a_->fixed() ) {
238  a_->fixed( kFALSE );
239  this->addFloatingParameter( a_ );
240  } else {
241  std::cerr << "WARNING in LauLASSRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
242  << std::endl;
243  }
244  } else if ( name == "r" ) {
245  if ( r_->fixed() ) {
246  r_->fixed( kFALSE );
247  this->addFloatingParameter( r_ );
248  } else {
249  std::cerr << "WARNING in LauLASSRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
250  << std::endl;
251  }
252  } else {
253  std::cerr << "WARNING in LauLASSRes::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
254  << std::endl;
255  }
256 }
257 
259 {
260  if ( name == "a" ) {
261  return a_;
262  } else if ( name == "r" ) {
263  return r_;
264  } else {
265  std::cerr << "WARNING in LauLASSRes::getResonanceParameter: Parameter name not reconised."
266  << std::endl;
267  return 0;
268  }
269 }
270 
271 void LauLASSRes::setEffectiveRange( const Double_t r )
272 {
273  r_->value( r );
274  r_->genValue( r );
275  r_->initValue( r );
276 }
277 
278 void LauLASSRes::setScatteringLength( const Double_t a )
279 {
280  a_->value( a );
281  a_->genValue( a );
282  a_->initValue( a );
283 }
File containing declaration of LauResonanceInfo class.
const TString & getSanitisedName() const
Get the name of the resonance.
virtual ~LauLASSRes()
Destructor.
Definition: LauLASSRes.cc:76
Double_t getWidth() const
Get the width of the resonance.
Double_t getScatteringLength() const
Get the scattering length range parameter.
Definition: LauLASSRes.hh:131
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
LauParameter * getMassPar()
Get the mass parameter of the resonance.
Double_t value() const
The value of the parameter.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
Double_t resMass_
The resonance mass.
Definition: LauLASSRes.hh:173
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
Definition: LauLASSRes.cc:116
virtual void initialise()
Initialise the model.
Definition: LauLASSRes.cc:80
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
Bool_t fixEffectiveRange() const
See if the effective range parameter is fixed or floating.
Definition: LauLASSRes.hh:137
void setEffectiveRange(const Double_t r)
Set the effective range parameter value.
Definition: LauLASSRes.cc:271
Int_t getSpin() const
Get the spin of the resonance.
Double_t resMassSq_
Square of the resonance mass.
Definition: LauLASSRes.hh:175
Double_t mDaugDiff_
Difference between the daughter masses.
Definition: LauLASSRes.hh:169
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
Definition: LauLASSRes.cc:234
Class for defining a complex number.
Definition: LauComplex.hh:61
void changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin)
Allow the mass, width and spin of the resonance to be changed.
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
Bool_t fixScatteringLength() const
See if the scattering length parameter is fixed or floating.
Definition: LauLASSRes.hh:143
Double_t mDaugDiffSq_
Square of mDaugDiff.
Definition: LauLASSRes.hh:171
void rescale(Double_t scaleVal)
Scale this by a factor.
Definition: LauComplex.hh:301
LauParameter * r_
LASS effective range parameter.
Definition: LauLASSRes.hh:178
File containing declaration of LauLASSRes class.
void setScatteringLength(const Double_t a)
Set the scattering length parameter value.
Definition: LauLASSRes.cc:278
Double_t cutOff_
LASS cut off.
Definition: LauLASSRes.hh:183
Bool_t fixWidth() const
Get the status of resonance width (fixed or released)
Double_t getMass() const
Get the mass of the resonance.
Double_t mDaugSum_
Sum of the daughter masses.
Definition: LauLASSRes.hh:165
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of a resonance parameter.
Definition: LauLASSRes.cc:217
void addExtraParameter(LauParameter *param, const Bool_t independentPar=kFALSE)
Add an extra parameter of the resonance.
LauParameter * a_
LASS scattering length parameter.
Definition: LauLASSRes.hh:180
Bool_t fixed() const
Check whether the parameter is fixed or floated.
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
Class for defining the properties of a resonant particle.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
Double_t getMassDaug1() const
Get the mass of daughter 1.
Double_t getQ() const
Get the current value of the daughter momentum in the resonance rest frame.
Bool_t fixMass() const
Get the status of resonance mass (fixed or released)
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Definition: LauLASSRes.cc:194
Double_t getMassDaug2() const
Get the mass of daughter 2.
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
Double_t getEffectiveRange() const
Get the effective range parameter.
Definition: LauLASSRes.hh:125
Double_t mDaugSumSq_
Square of mDaugSum.
Definition: LauLASSRes.hh:167
Double_t genValue() const
The value generated for the parameter.
Double_t q0_
Decay momentum of either daughter in the resonance rest frame.
Definition: LauLASSRes.hh:163
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Definition: LauLASSRes.cc:258
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
void calcQ0()
Utility function to calculate the q0 value.
Definition: LauLASSRes.cc:104
LauParameter * getWidthPar()
Get the width parameter of the resonance.
Double_t initValue() const
The initial value of the parameter.
LauLASSRes(LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
Definition: LauLASSRes.cc:36