laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauSigmaRes.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 "LauSigmaRes.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  mPiSq4_( 4.0 * LauConstants::mPiSq ),
41  sAdler_( LauConstants::mPiSq * 0.5 ),
42  b1_( 0 ),
43  b2_( 0 ),
44  a_( 0 ),
45  m0_( 0 )
46 {
47  // Initialise various constants
49  sAdler_ = LauConstants::mPiSq * 0.5; // Adler zero at 0.5*(mpi)^2
50 
51  // constant factors from BES data
52  const Double_t b1Val = 0.5843;
53  const Double_t b2Val = 1.6663;
54  const Double_t aVal = 1.082;
55  const Double_t m0Val = 0.9264;
56 
57  const TString& parNameBase = this->getSanitisedName();
58 
59  TString b1Name( parNameBase );
60  b1Name += "_b1";
61  b1_ = resInfo->getExtraParameter( b1Name );
62  if ( b1_ == 0 ) {
63  b1_ = new LauParameter( b1Name, b1Val, 0.0, 100.0, kTRUE );
64  b1_->secondStage( kTRUE );
65  resInfo->addExtraParameter( b1_ );
66  }
67 
68  TString b2Name( parNameBase );
69  b2Name += "_b2";
70  b2_ = resInfo->getExtraParameter( b2Name );
71  if ( b2_ == 0 ) {
72  b2_ = new LauParameter( b2Name, b2Val, 0.0, 100.0, kTRUE );
73  b2_->secondStage( kTRUE );
74  resInfo->addExtraParameter( b2_ );
75  }
76 
77  TString aName( parNameBase );
78  aName += "_A";
79  a_ = resInfo->getExtraParameter( aName );
80  if ( a_ == 0 ) {
81  a_ = new LauParameter( aName, aVal, 0.0, 10.0, kTRUE );
82  a_->secondStage( kTRUE );
83  resInfo->addExtraParameter( a_ );
84  }
85 
86  TString m0Name( parNameBase );
87  m0Name += "_m0";
88  m0_ = resInfo->getExtraParameter( m0Name );
89  if ( m0_ == 0 ) {
90  m0_ = new LauParameter( m0Name, m0Val, 0.0, 10.0, kTRUE );
91  m0_->secondStage( kTRUE );
92  resInfo->addExtraParameter( m0_ );
93  }
94 }
95 
97 {
98 }
99 
101 {
102  this->checkDaughterTypes();
103 
104  Double_t resSpin = this->getSpin();
105  if ( resSpin != 0 ) {
106  std::cerr << "WARNING in LauSigmaRes::initialise : Resonance spin is " << resSpin << "."
107  << std::endl;
108  std::cerr << " : Sigma amplitude is only for scalers, resetting spin to 0."
109  << std::endl;
110  this->changeResonance( -1.0, -1.0, 0 );
111  }
112 }
113 
115 {
116  // Check that the daughter tracks are the same type. Otherwise issue a warning
117  // and set the type to be pion for the Sigma distribution.
118  Int_t resPairAmpInt = this->getPairInt();
119  if ( resPairAmpInt < 1 || resPairAmpInt > 3 ) {
120  std::cerr << "WARNING in LauSigmaRes::checkDaughterTypes : resPairAmpInt = " << resPairAmpInt
121  << " is out of the range [1,2,3]." << std::endl;
122  return;
123  }
124 
125  const TString& nameDaug1 = this->getNameDaug1();
126  const TString& nameDaug2 = this->getNameDaug2();
127  if ( ! nameDaug1.CompareTo( nameDaug2, TString::kExact ) ) {
128  // Daughter types agree. Find out if we have pion or kaon system
129  if ( ! nameDaug1.Contains( "pi" ) ) {
130  std::cerr << "ERROR in LauSigmaRes::checkDaughterTypes : Sigma model is using daughters \""
131  << nameDaug1 << "\" and \"" << nameDaug2 << "\", which are not pions."
132  << std::endl;
133  }
134  }
135 }
136 
137 LauComplex LauSigmaRes::resAmp( Double_t mass, Double_t spinTerm )
138 {
139  // This function returns the complex dynamical amplitude for a Sigma distribution
140  // given the invariant mass and cos(helicity) values.
141 
142  // First check that the appropriate daughters are either pi+pi- or K+K-
143  // Check that the daughter tracks are the same type. Otherwise issue a warning
144  // and set the type to be pion for the Sigma distribution. Returns the
145  // integer defined by the enum LauSigmaRes::SigmaPartType.
146 
147  Double_t s = mass * mass; // Invariant mass squared combination for the system
148  Double_t rho( 0.0 ); // Phase-space factor
149  if ( s > mPiSq4_ ) {
150  rho = TMath::Sqrt( 1.0 - mPiSq4_ / s );
151  }
152 
153  const Double_t m0Val = this->getM0Value();
154  const Double_t m0Sq = m0Val * m0Val;
155 
156  const Double_t aVal = this->getAValue();
157  const Double_t b1Val = this->getB1Value();
158  const Double_t b2Val = this->getB2Value();
159 
160  Double_t f = b2Val * s + b1Val; // f(s) function
161  Double_t numerator = s - sAdler_;
162  Double_t denom = m0Sq - sAdler_;
163  Double_t gamma( 0.0 );
164  if ( TMath::Abs( denom ) > 1e-10 && TMath::Abs( aVal ) > 1e-10 ) {
165  // Decay width of the system
166  gamma = rho * ( numerator / denom ) * f * TMath::Exp( -( s - m0Sq ) / aVal );
167  }
168 
169  // Now form the complex amplitude - use relativistic BW form (without barrier factors)
170  // Note that the M factor in the denominator is not the "pole" at ~500 MeV, but is
171  // m0_ = 0.9264, the mass when the phase shift goes through 90 degrees.
172 
173  Double_t dMSq = m0Sq - s;
174  Double_t widthTerm = gamma * m0Val;
175  LauComplex resAmplitude( dMSq, widthTerm );
176 
177  Double_t denomFactor = dMSq * dMSq + widthTerm * widthTerm;
178 
179  Double_t invDenomFactor = 0.0;
180  if ( denomFactor > 1e-10 ) {
181  invDenomFactor = 1.0 / denomFactor;
182  }
183 
184  resAmplitude.rescale( spinTerm * invDenomFactor );
185 
186  return resAmplitude;
187 }
188 
189 const std::vector<LauParameter*>& LauSigmaRes::getFloatingParameters()
190 {
191  this->clearFloatingParameters();
192 
193  if ( ! this->fixB1Value() ) {
194  this->addFloatingParameter( b1_ );
195  }
196 
197  if ( ! this->fixB2Value() ) {
198  this->addFloatingParameter( b2_ );
199  }
200 
201  if ( ! this->fixAValue() ) {
202  this->addFloatingParameter( a_ );
203  }
204 
205  if ( ! this->fixM0Value() ) {
206  this->addFloatingParameter( m0_ );
207  }
208 
209  return this->getParameters();
210 }
211 
212 void LauSigmaRes::setResonanceParameter( const TString& name, const Double_t value )
213 {
214  // Set various parameters for the lineshape
215  if ( name == "b1" ) {
216  this->setB1Value( value );
217  std::cout << "INFO in LauSigmaRes::setResonanceParameter : Setting parameter b1 = "
218  << this->getB1Value() << std::endl;
219  } else if ( name == "b2" ) {
220  this->setB2Value( value );
221  std::cout << "INFO in LauSigmaRes::setResonanceParameter : Setting parameter b2 = "
222  << this->getB2Value() << std::endl;
223  } else if ( name == "A" ) {
224  this->setAValue( value );
225  std::cout << "INFO in LauSigmaRes::setResonanceParameter : Setting parameter A = "
226  << this->getAValue() << std::endl;
227  } else if ( name == "m0" ) {
228  this->setM0Value( value );
229  std::cout << "INFO in LauSigmaRes::setResonanceParameter : Setting parameter m0 = "
230  << this->getM0Value() << std::endl;
231  } else {
232  std::cerr << "WARNING in LauSigmaRes::setResonanceParameter: Parameter name not reconised. No parameter changes made."
233  << std::endl;
234  }
235 }
236 
238 {
239  if ( name == "b1" ) {
240  if ( b1_->fixed() ) {
241  b1_->fixed( kFALSE );
242  this->addFloatingParameter( b1_ );
243  } else {
244  std::cerr << "WARNING in LauSigmaRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
245  << std::endl;
246  }
247  } else if ( name == "b2" ) {
248  if ( b2_->fixed() ) {
249  b2_->fixed( kFALSE );
250  this->addFloatingParameter( b2_ );
251  } else {
252  std::cerr << "WARNING in LauSigmaRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
253  << std::endl;
254  }
255  } else if ( name == "A" ) {
256  if ( a_->fixed() ) {
257  a_->fixed( kFALSE );
258  this->addFloatingParameter( a_ );
259  } else {
260  std::cerr << "WARNING in LauSigmaRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
261  << std::endl;
262  }
263  } else if ( name == "m0" ) {
264  if ( m0_->fixed() ) {
265  m0_->fixed( kFALSE );
266  this->addFloatingParameter( m0_ );
267  } else {
268  std::cerr << "WARNING in LauSigmaRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
269  << std::endl;
270  }
271  } else {
272  std::cerr << "WARNING in LauSigmaRes::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
273  << std::endl;
274  }
275 }
276 
278 {
279  if ( name == "b1" ) {
280  return b1_;
281  } else if ( name == "b2" ) {
282  return b2_;
283  } else if ( name == "A" ) {
284  return a_;
285  } else if ( name == "m0" ) {
286  return m0_;
287  } else {
288  std::cerr << "WARNING in LauSigmaRes::getResonanceParameter: Parameter name not reconised."
289  << std::endl;
290  return 0;
291  }
292 }
293 
294 void LauSigmaRes::setB1Value( const Double_t b1 )
295 {
296  b1_->value( b1 );
297  b1_->genValue( b1 );
298  b1_->initValue( b1 );
299 }
300 
301 void LauSigmaRes::setB2Value( const Double_t b2 )
302 {
303  b2_->value( b2 );
304  b2_->genValue( b2 );
305  b2_->initValue( b2 );
306 }
307 
308 void LauSigmaRes::setAValue( const Double_t A )
309 {
310  a_->value( A );
311  a_->genValue( A );
312  a_->initValue( A );
313 }
314 
315 void LauSigmaRes::setM0Value( const Double_t m0 )
316 {
317  m0_->value( m0 );
318  m0_->genValue( m0 );
319  m0_->initValue( m0 );
320 }
Double_t sAdler_
Defined as 0.5*mPi*mPi.
Definition: LauSigmaRes.hh:189
File containing declaration of LauSigmaRes class.
File containing declaration of LauResonanceInfo class.
const Double_t mPiSq
Square of pi+- mass.
Definition: LauConstants.hh:80
const TString & getSanitisedName() const
Get the name of the resonance.
LauSigmaRes(LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
Definition: LauSigmaRes.cc:36
virtual void initialise()
Initialise the model.
Definition: LauSigmaRes.cc:100
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Definition: LauSigmaRes.cc:212
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: LauSigmaRes.cc:277
void setB2Value(const Double_t b2)
Set the b2 parameter.
Definition: LauSigmaRes.cc:301
void setM0Value(const Double_t m0)
Set the m0 parameter.
Definition: LauSigmaRes.cc:315
void setB1Value(const Double_t b1)
Set the b1 parameter.
Definition: LauSigmaRes.cc:294
LauParameter * m0_
Factor from BES data.
Definition: LauSigmaRes.hh:198
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
Double_t mPiSq4_
Defined as 4*mPi*mPi.
Definition: LauSigmaRes.hh:187
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
LauParameter * a_
Factor from BES data.
Definition: LauSigmaRes.hh:196
Int_t getSpin() const
Get the spin of the resonance.
Bool_t fixB1Value() const
Fix the b1 parameter value.
Definition: LauSigmaRes.hh:149
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit.
void setAValue(const Double_t A)
Set the A parameter.
Definition: LauSigmaRes.cc:308
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.
Namespace to contain various constants.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
virtual ~LauSigmaRes()
Destructor.
Definition: LauSigmaRes.cc:96
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant ampltiude.
Definition: LauSigmaRes.cc:137
void rescale(Double_t scaleVal)
Scale this by a factor.
Definition: LauComplex.hh:301
Double_t getB1Value() const
Get the b1 parameter value.
Definition: LauSigmaRes.hh:125
Double_t getAValue() const
Get the A parameter value.
Definition: LauSigmaRes.hh:137
LauParameter * b1_
Factor from BES data.
Definition: LauSigmaRes.hh:192
void checkDaughterTypes() const
Check that both daughters are the same type of particle.
Definition: LauSigmaRes.cc:114
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.
LauParameter()
Default constructor.
Definition: LauParameter.cc:40
Bool_t fixM0Value() const
Fix the m0 parameter value.
Definition: LauSigmaRes.hh:167
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Definition: LauSigmaRes.cc:189
Class for defining the properties of a resonant particle.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
Double_t getB2Value() const
Get the b2 parameter value.
Definition: LauSigmaRes.hh:131
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
Definition: LauSigmaRes.cc:237
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
Bool_t fixB2Value() const
Fix the b2 parameter value.
Definition: LauSigmaRes.hh:155
Bool_t fixAValue() const
Fix the A parameter value.
Definition: LauSigmaRes.hh:161
Double_t genValue() const
The value generated for the parameter.
LauParameter * b2_
Factor from BES data.
Definition: LauSigmaRes.hh:194
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
TString getNameDaug1() const
Get the name of the first daughter of the resonance.
Double_t getM0Value() const
Get the m0 parameter value.
Definition: LauSigmaRes.hh:143
TString getNameDaug2() const
Get the name of the second daughter of the resonance.
Double_t initValue() const
The initial value of the parameter.