laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauRescatteringRes.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2018 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 "LauRescatteringRes.hh"
30 
31 #include "LauConstants.hh"
32 #include "LauDaughters.hh"
33 #include "LauParameter.hh"
34 #include "LauResonanceInfo.hh"
35 
36 #include "TMath.h"
37 #include "TSystem.h"
38 
39 #include <iostream>
40 
43  const Int_t resPairAmpInt,
44  const LauDaughters* daughters ) :
45  LauAbsResonance( resInfo, resPairAmpInt, daughters ),
46  lambdaPiPi_( 0 ),
47  lambdaKK_( 0 ),
48  Mf_( 0 ),
49  Ms_( 0 ),
50  Mprime_( 0 ),
51  Eps1_( 0 ),
52  Eps2_( 0 ),
53  C0_( 0 ),
54  model_( resType )
55 {
56  TString parNameBase = this->getSanitisedName();
57  TString lambdaPiPiName( parNameBase );
58  lambdaPiPiName += "_lambdaPiPi";
59  lambdaPiPi_ = resInfo->getExtraParameter( lambdaPiPiName );
60  if ( lambdaPiPi_ == 0 ) {
61  lambdaPiPi_ = new LauParameter( lambdaPiPiName, 1.0, 0.0, 10.0, kTRUE );
62  lambdaPiPi_->secondStage( kTRUE );
63  resInfo->addExtraParameter( lambdaPiPi_ );
64  }
65  TString lambdaKKName( parNameBase );
66  lambdaKKName += "_lambdaKK";
67  lambdaKK_ = resInfo->getExtraParameter( lambdaKKName );
68  if ( lambdaKK_ == 0 ) {
69  lambdaKK_ = new LauParameter( lambdaKKName, 2.8, 0.0, 10.0, kTRUE );
70  lambdaKK_->secondStage( kTRUE );
71  resInfo->addExtraParameter( lambdaKK_ );
72  }
73  TString MfName( parNameBase );
74  MfName += "_Mf";
75  Mf_ = resInfo->getExtraParameter( MfName );
76  if ( Mf_ == 0 ) {
77  Mf_ = new LauParameter( MfName, 1.32, 0.0, 10.0, kTRUE );
78  Mf_->secondStage( kTRUE );
79  resInfo->addExtraParameter( Mf_ );
80  }
81  TString MsName( parNameBase );
82  MsName += "_Ms";
83  Ms_ = resInfo->getExtraParameter( MsName );
84  if ( Ms_ == 0 ) {
85  Ms_ = new LauParameter( MsName, 0.92, 0.0, 10.0, kTRUE );
86  Ms_->secondStage( kTRUE );
87  resInfo->addExtraParameter( Ms_ );
88  }
89  TString MprimeName( parNameBase );
90  MprimeName += "_Mprime";
91  Mprime_ = resInfo->getExtraParameter( MprimeName );
92  if ( Mprime_ == 0 ) {
93  Mprime_ = new LauParameter( MprimeName, 1.5, 0.0, 10.0, kTRUE );
94  Mprime_->secondStage( kTRUE );
95  resInfo->addExtraParameter( Mprime_ );
96  }
97 
98  TString Eps1Name( parNameBase );
99  Eps1Name += "_Eps1";
100  Eps1_ = resInfo->getExtraParameter( Eps1Name );
101  if ( Eps1_ == 0 ) {
102  Eps1_ = new LauParameter( Eps1Name, 2.4, 0.0, 10.0, kTRUE );
103  Eps1_->secondStage( kTRUE );
104  resInfo->addExtraParameter( Eps1_ );
105  }
106 
107  TString Eps2Name( parNameBase );
108  Eps2Name += "_Eps2";
109  Eps2_ = resInfo->getExtraParameter( Eps2Name );
110  if ( Eps2_ == 0 ) {
111  Eps2_ = new LauParameter( Eps2Name, -5.5, -10.0, 10.0, kTRUE );
112  Eps2_->secondStage( kTRUE );
113  resInfo->addExtraParameter( Eps2_ );
114  }
115 
116  TString C0Name( parNameBase );
117  C0Name += "_C0";
118  C0_ = resInfo->getExtraParameter( C0Name );
119  if ( C0_ == 0 ) {
120  C0_ = new LauParameter( C0Name, 1.3, 0.0, 10.0, kTRUE );
121  C0_->secondStage( kTRUE );
122  resInfo->addExtraParameter( C0_ );
123  }
124 }
125 
127 {
128 }
129 
131 {
132  const LauDaughters* daughters = this->getDaughters();
133  Int_t resPairAmpInt = this->getPairInt();
134  if ( daughters->gotSymmetricalDP() && resPairAmpInt != 3 ) {
135  std::cerr << "WARNING in LauRescatteringRes::initialise : Dalitz plot is symmetric - this lineshape is not appropriate."
136  << std::endl;
137  std::cerr << "WARNING I think that this warning is not appropiate because Laura Simetrize at LauIsobarModel level."
138  << std::endl;
139  }
140 
143  std::cerr << "WARNING in LauRescatteringRes::initialise : Unknown model requested, defaulting to exponential."
144  << std::endl;
146  }
147 
148  if ( ( model_ != LauAbsResonance::RescatteringNoInter ) && ( this->getSpin() != 0 ) ) {
149  std::cerr << "WARNING in LauRescatteringRes::initialise : Non-zero spin will be ignored for this model - perhaps you should use LauAbsResonance::BelleSymNRNoInter instead"
150  << std::endl;
151  }
152 }
153 
155 {
156  // This function returns the complex dynamical amplitude for a Reescatering distribution o original
157  // pelaez paper parameters Eq 2.15a [Pelaez et Yndúrain: arXiv:hep-ph/0411334v2 Mar 2005]
158  Double_t Mprime = this->getMprime();
159  Double_t Mf = this->getMf();
160  Double_t Ms = this->getMs();
161  Double_t eps1 = this->getEps1();
162  Double_t eps2 = this->getEps2();
163  Double_t c0 = this->getC0();
164  Double_t lambPiPi = this->getLambdaPiPi();
165  Double_t lambKK = this->getLambdaKK();
166 
167  Double_t mk = LauConstants::mK;
168 
169  // Calculate Mandelstam variables: s = m_13^2, t = m_23^2, u = m_12^2.
170  Double_t s = 0;
171  Double_t t = 0;
172  Int_t resPairAmpInt = getPairInt();
173 
174  if ( resPairAmpInt == 1 ) {
175  s = kinematics->getm23Sq();
176  t = kinematics->getm12Sq();
177  } else if ( resPairAmpInt == 2 ) {
178  s = kinematics->getm13Sq();
179  t = kinematics->getm23Sq();
180  } else if ( resPairAmpInt == 3 ) {
181  s = kinematics->getm12Sq();
182  t = kinematics->getm13Sq();
183  } else {
184  std::cerr << "ERROR in LauAbsResonance::amplitude : Nonsense setup of resPairAmp array."
185  << std::endl;
186  gSystem->Exit( EXIT_FAILURE );
187  }
188 
189  // Calculate amplitude for s variable.
190  Double_t mass_s = TMath::Sqrt( s );
191  Double_t k2Square_s = ( s - 4.0 * mk * mk ) / 4.0;
192  Double_t k2Abs_s = 0;
193  if ( k2Square_s > 0 )
194  k2Abs_s = TMath::Sqrt( k2Square_s );
195  else
196  k2Abs_s = TMath::Sqrt( -1.0 * k2Square_s );
197  Double_t cotdelta0_s = c0 * ( s - Ms * Ms ) * ( Mf * Mf - s ) * k2Abs_s /
198  ( Mf * Mf * mass_s * k2Square_s ); // Eq 2.15a
199  Double_t delta0_s = TMath::ATan( 1.0 / cotdelta0_s );
200  Double_t eta0_s = 1.0 - ( eps1 * k2Abs_s / mass_s + eps2 * k2Square_s / s ) *
201  ( Mprime * Mprime - s ) / s; // Eq 2.15a
202 
203  if ( ( mass_s < 2.0 * mk ) || ( mass_s > Mprime ) )
204  eta0_s = 1;
205  Double_t mag_s = TMath::Sqrt( 1 - eta0_s * eta0_s );
206 
207  Double_t tauRe_s = mag_s * TMath::Cos( 2.0 * delta0_s );
208  Double_t tauIm_s = mag_s * TMath::Sin( 2.0 * delta0_s );
209 
210  Double_t NR1_s = 1.0 / ( 1.0 + s / ( lambPiPi * lambPiPi ) );
211  Double_t NR2_s = 1.0 / ( 1.0 + s / ( lambKK * lambKK ) );
212 
213  //LauComplex resAmplitude(-tauIm_s*NR1_s*NR2_s - tauIm_t*NR1_t*NR2_t, tauRe_s*NR1_s*NR2_s + tauRe_t*NR1_t*NR2_t );
214  if ( ( model_ == LauAbsResonance::RescatteringNoInter ) && ( t <= s ) ) {
215  NR1_s = 0.0;
216  NR2_s = 0.0;
217  tauRe_s = 0.0;
218  tauIm_s = 0.0;
219  }
220 
221  LauComplex resAmplitude( -tauIm_s * NR1_s * NR2_s, tauRe_s * NR1_s * NR2_s );
222 
223  return resAmplitude;
224 }
225 
226 LauComplex LauRescatteringRes::resAmp( Double_t mass, Double_t spinTerm )
227 {
228  std::cerr << "ERROR in LauRescatteringRes : This method should never be called." << std::endl;
229  std::cerr << " : Returning zero amplitude for mass = " << mass
230  << " and spinTerm = " << spinTerm << "." << std::endl;
231  return LauComplex( 0.0, 0.0 );
232 }
233 
234 const std::vector<LauParameter*>& LauRescatteringRes::getFloatingParameters()
235 {
236  this->clearFloatingParameters();
237 
238  if ( ! this->fixLambdaPiPi() ) {
240  }
241  if ( ! this->fixLambdaKK() ) {
243  }
244  if ( ! this->fixMf() ) {
245  this->addFloatingParameter( Mf_ );
246  }
247  if ( ! this->fixMs() ) {
248  this->addFloatingParameter( Ms_ );
249  }
250  if ( ! this->fixMprime() ) {
251  this->addFloatingParameter( Mprime_ );
252  }
253  if ( ! this->fixEps1() ) {
254  this->addFloatingParameter( Eps1_ );
255  }
256  if ( ! this->fixEps2() ) {
257  this->addFloatingParameter( Eps2_ );
258  }
259  if ( ! this->fixC0() ) {
260  this->addFloatingParameter( C0_ );
261  }
262 
263  return this->getParameters();
264 }
265 
266 void LauRescatteringRes::setResonanceParameter( const TString& name, const Double_t value )
267 {
268  // Set various parameters for the lineshape
269  if ( name == "lambdaPiPi" ) {
270  this->setLambdaPiPi( value );
271  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter lambdaPiPi = "
272  << this->getLambdaPiPi() << std::endl;
273  } else if ( name == "lambdaKK" ) {
274  this->setLambdaKK( value );
275  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter lambdaKK = "
276  << this->getLambdaKK() << std::endl;
277  } else if ( name == "Mf" ) {
278  this->setMf( value );
279  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter Mf = "
280  << this->getMf() << std::endl;
281  } else if ( name == "Ms" ) {
282  this->setMs( value );
283  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter Ms = "
284  << this->getMs() << std::endl;
285  } else if ( name == "Mprime" ) {
286  this->setMprime( value );
287  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter Mprime = "
288  << this->getMprime() << std::endl;
289  } else if ( name == "Eps1" ) {
290  this->setEps1( value );
291  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter eps1 = "
292  << this->getEps1() << std::endl;
293  } else if ( name == "Eps2" ) {
294  this->setEps2( value );
295  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter eps2 = "
296  << this->getEps2() << std::endl;
297  } else if ( name == "C0" ) {
298  this->setC0( value );
299  std::cout << "INFO in LauRescatteringRes::setResonanceParameter : Setting parameter eps2 = "
300  << this->getC0() << std::endl;
301  } else {
302  std::cerr << "WARNING in LauRescatteringRes::setResonanceParameter: Parameter name not reconised. No parameter changes made."
303  << std::endl;
304  }
305 }
306 
308 {
309  if ( name == "lambdaPiPi" ) {
310  if ( lambdaPiPi_->fixed() ) {
311  lambdaPiPi_->fixed( kFALSE );
313  } else {
314  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
315  << std::endl;
316  }
317  } else if ( name == "lambdaKK" ) {
318  if ( lambdaKK_->fixed() ) {
319  lambdaKK_->fixed( kFALSE );
321  } else {
322  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
323  << std::endl;
324  }
325  } else if ( name == "Mf" ) {
326  if ( Mf_->fixed() ) {
327  Mf_->fixed( kFALSE );
328  this->addFloatingParameter( Mf_ );
329  } else {
330  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
331  << std::endl;
332  }
333  } else if ( name == "Ms" ) {
334  if ( Ms_->fixed() ) {
335  Ms_->fixed( kFALSE );
336  this->addFloatingParameter( Ms_ );
337  } else {
338  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
339  << std::endl;
340  }
341  } else if ( name == "Mprime" ) {
342  if ( Mprime_->fixed() ) {
343  Mprime_->fixed( kFALSE );
344  this->addFloatingParameter( Mprime_ );
345  } else {
346  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
347  << std::endl;
348  }
349  } else if ( name == "Eps1" ) {
350  if ( Eps1_->fixed() ) {
351  Eps1_->fixed( kFALSE );
352  this->addFloatingParameter( Eps1_ );
353  } else {
354  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
355  << std::endl;
356  }
357  } else if ( name == "Eps2" ) {
358  if ( Eps2_->fixed() ) {
359  Eps2_->fixed( kFALSE );
360  this->addFloatingParameter( Eps2_ );
361  } else {
362  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
363  << std::endl;
364  }
365  } else if ( name == "C0" ) {
366  if ( C0_->fixed() ) {
367  C0_->fixed( kFALSE );
368  this->addFloatingParameter( C0_ );
369  } else {
370  std::cerr << "WARNING in LauRescatteringRes::floatResonanceParameter: Parameter already floating. No parameter changes made."
371  << std::endl;
372  }
373  }
374 
375  else {
376  std::cerr << "WARNING in LauRescatteringRes::fixResonanceParameter: Parameter name not reconised. No parameter changes made."
377  << std::endl;
378  }
379 }
380 
382 {
383  if ( name == "lambdaPiPi" ) {
384  return lambdaPiPi_;
385  } else if ( name == "lambdaKK" ) {
386  return lambdaKK_;
387  }
388  if ( name == "Mf" ) {
389  return Mf_;
390  } else if ( name == "Ms" ) {
391  return Ms_;
392  }
393  if ( name == "Mprime" ) {
394  return Mprime_;
395  }
396  if ( name == "Eps1" ) {
397  return Eps1_;
398  }
399  if ( name == "Eps2" ) {
400  return Eps2_;
401  }
402  if ( name == "C0" ) {
403  return C0_;
404  } else {
405  std::cerr << "WARNING in LauRescatteringRes::getResonanceParameter: Parameter name not reconised."
406  << std::endl;
407  return 0;
408  }
409 }
410 
411 void LauRescatteringRes::setLambdaPiPi( const Double_t lambda )
412 {
413  lambdaPiPi_->value( lambda );
414  lambdaPiPi_->genValue( lambda );
415  lambdaPiPi_->initValue( lambda );
416 }
417 
418 void LauRescatteringRes::setLambdaKK( const Double_t lambda )
419 {
420  lambdaKK_->value( lambda );
421  lambdaKK_->genValue( lambda );
422  lambdaKK_->initValue( lambda );
423 }
424 
425 void LauRescatteringRes::setMf( const Double_t Mf )
426 {
427  Mf_->value( Mf );
428  Mf_->genValue( Mf );
429  Mf_->initValue( Mf );
430 }
431 
432 void LauRescatteringRes::setMs( const Double_t Ms )
433 {
434  Ms_->value( Ms );
435  Ms_->genValue( Ms );
436  Ms_->initValue( Ms );
437 }
438 
439 void LauRescatteringRes::setMprime( const Double_t Mprime )
440 {
441  Mprime_->value( Mprime );
442  Mprime_->genValue( Mprime );
443  Mprime_->initValue( Mprime );
444 }
445 
446 void LauRescatteringRes::setEps1( const Double_t Eps1 )
447 {
448  Eps1_->value( Eps1 );
449  Eps1_->genValue( Eps1 );
450  Eps1_->initValue( Eps1 );
451 }
452 
453 void LauRescatteringRes::setEps2( const Double_t Eps2 )
454 {
455  Eps2_->value( Eps2 );
456  Eps2_->genValue( Eps2 );
457  Eps2_->initValue( Eps2 );
458 }
459 
460 void LauRescatteringRes::setC0( const Double_t C0 )
461 {
462  C0_->value( C0 );
463  C0_->genValue( C0 );
464  C0_->initValue( C0 );
465 }
virtual LauComplex amplitude(const LauKinematics *kinematics)
Get the complex dynamical amplitude.
File containing declaration of LauResonanceInfo class.
LauAbsResonance::LauResonanceModel model_
The model to use.
const TString & getSanitisedName() const
Get the name of the resonance.
LauParameter * C0_
the term for the C0
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
void setLambdaKK(const Double_t lambda)
Set the parameter lambdaKK, the term for the KK.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
File containing declaration of LauRescatteringRes class.
Double_t value() const
The value of the parameter.
virtual const std::vector< LauParameter * > & getFloatingParameters()
Retrieve the resonance parameters, e.g. so that they can be loaded into a fit.
Double_t getm12Sq() const
Get the m12 invariant mass square.
Bool_t fixMprime() const
See if the Mprime parameter is fixed or floating.
Bool_t fixMs() const
See if the Ms parameter is fixed or floating.
std::vector< LauParameter * > & getParameters()
Access the list of floating parameters.
File containing declaration of LauParameter class.
Double_t getC0() const
Get the C0.
LauParameter * Eps2_
the term for the Eps2
void setEps1(const Double_t Eps1)
Set the parameter Eps1.
LauParameter * Mf_
the term for the Mf_
LauParameter * getExtraParameter(const TString &parName)
Retrieve an extra parameter of the resonance.
LauRescatteringRes(LauResonanceInfo *resInfo, const LauAbsResonance::LauResonanceModel resType, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor.
void setMprime(const Double_t Mprime)
Set the parameter Mprime.
Double_t getLambdaPiPi() const
Get the lambdaPiPi, the term for the PiPi.
Int_t getSpin() const
Get the spin of the resonance.
void setMf(const Double_t Mf)
Set the parameter Mf.
Bool_t secondStage() const
Check whether the parameter should be floated only in the second stage of a two stage fit.
Bool_t fixEps2() const
See if the Eps2 parameter is fixed or floating.
Double_t getMf() const
Get the Mf.
Class for defining a complex number.
Definition: LauComplex.hh:61
Double_t getMprime() const
Get the Mprime.
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
Bool_t fixLambdaPiPi() const
See if the lambdaPiPi parameter is fixed or floating.
Int_t getPairInt() const
Get the integer to identify which DP axis the resonance belongs to.
const LauDaughters * getDaughters() const
Access the daughters object.
Bool_t fixLambdaKK() const
See if the lambdaKK parameter is fixed or floating.
Double_t getMs() const
Get the Ms.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
virtual void initialise()
Initialise the model.
Double_t getEps2() const
Get the Eps2.
Bool_t fixEps1() const
See if the Eps1 parameter is fixed or floating.
File containing declaration of LauDaughters class.
LauParameter * Eps1_
the term for the Eps1
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
Class for defining the properties of a resonant particle.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
LauParameter * lambdaPiPi_
the term for the PiPi
LauParameter * lambdaKK_
the term for the KK
void setEps2(const Double_t Eps2)
Set the parameter Eps2.
const Double_t mK
Mass of K+- (GeV/c^2)
Definition: LauConstants.hh:59
LauParameter * Ms_
the term for the Ms
LauParameter * Mprime_
the term for the Mprime
void setLambdaPiPi(const Double_t lambda)
Set the parameter lambdaPiPi, the term for the PiPi.
Double_t getLambdaKK() const
Get the lambdaKK, the term for the KK.
Double_t getEps1() const
Get the Eps1.
void clearFloatingParameters()
Clear list of floating parameters.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc....
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)
Complex resonant amplitude.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Class for calculating 3-body kinematic quantities.
Bool_t fixC0() const
See if the C0 parameter is fixed or floating.
Double_t genValue() const
The value generated for the parameter.
virtual ~LauRescatteringRes()
Destructor.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getm23Sq() const
Get the m23 invariant mass square.
Double_t getm13Sq() const
Get the m13 invariant mass square.
Bool_t gotSymmetricalDP() const
Is Dalitz plot symmetric, i.e. 2 identical particles.
Definition: LauDaughters.hh:84
void setMs(const Double_t Ms)
Set the parameter Ms.
void setC0(const Double_t C0)
Set the parameter C0.
Double_t initValue() const
The initial value of the parameter.
LauResonanceModel
Define the allowed resonance types.
Bool_t fixMf() const
See if the Mf parameter is fixed or floating.