laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauAbsResonance.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 "LauAbsResonance.hh"
30 
31 #include "LauConstants.hh"
32 #include "LauDaughters.hh"
33 #include "LauKinematics.hh"
34 #include "LauParameter.hh"
35 #include "LauResonanceInfo.hh"
36 
37 #include "TSystem.h"
38 
39 #include <iostream>
40 
42 {
43  switch ( model ) {
44  case BW :
45  case RelBW :
46  case GS :
47  case Flatte :
48  case Sigma :
49  case Kappa :
50  case Dabba :
51  case LASS :
52  case LASS_BW :
53  case LASS_NR :
54  case EFKLLM :
55  case KMatrix :
56  case FlatNR :
57  case NRModel :
58  case BelleNR :
59  case PowerLawNR :
60  case BelleSymNR :
61  case BelleSymNRNoInter :
62  case TaylorNR :
63  case PolNR :
64  case Pole :
65  case PolarFFNR :
66  case PolarFFSymNR :
67  case PolarFFSymNRNoInter :
68  case Rescattering :
69  case Rescattering2 :
70  case RescatteringNoInter :
71  case MIPW_MagPhase :
72  case MIPW_RealImag :
73  case RhoOmegaMix_GS :
74  case RhoOmegaMix_RBW :
75  case RhoOmegaMix_GS_1 :
76  case RhoOmegaMix_RBW_1 :
77  break;
78  case GaussIncoh :
79  return true;
80  }
81  return false;
82 }
83 
84 // Constructor
86  const Int_t resPairAmpInt,
87  const LauDaughters* daughters ) :
88  resInfo_( resInfo ),
89  daughters_( daughters ),
90  resName_( ( resInfo != 0 ) ? resInfo->getName() : "" ),
91  sanitisedName_( ( resInfo != 0 ) ? resInfo->getSanitisedName() : "" ),
92  resMass_( ( resInfo != 0 ) ? resInfo->getMass() : 0 ),
93  resWidth_( ( resInfo != 0 ) ? resInfo->getWidth() : 0 ),
94  resSpin_( ( resInfo != 0 ) ? resInfo->getSpin() : 0 ),
95  resCharge_( ( resInfo != 0 ) ? resInfo->getCharge() : 0 ),
96  resPairAmpInt_( resPairAmpInt )
97 {
98  if ( resInfo == 0 ) {
99  std::cerr << "ERROR in LauAbsResonance constructor : null LauResonanceInfo object provided"
100  << std::endl;
101  gSystem->Exit( EXIT_FAILURE );
102  }
103 
104  if ( daughters_ == 0 ) {
105  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided"
106  << std::endl;
107  gSystem->Exit( EXIT_FAILURE );
108  }
109 
110  nameParent_ = this->getNameParent();
111  nameDaug1_ = this->getNameDaug1();
112  nameDaug2_ = this->getNameDaug2();
113  nameBachelor_ = this->getNameBachelor();
114  massParent_ = this->getMassParent();
115  massDaug1_ = this->getMassDaug1();
116  massDaug2_ = this->getMassDaug2();
117  massBachelor_ = this->getMassBachelor();
118  chargeParent_ = this->getChargeParent();
119  chargeDaug1_ = this->getChargeDaug1();
120  chargeDaug2_ = this->getChargeDaug2();
122 
123  // check that the total charge adds up to that of the resonance
124  Int_t totalCharge = chargeDaug1_ + chargeDaug2_;
125  if ( ( totalCharge != resCharge_ ) && ( resPairAmpInt_ != 0 ) ) {
126  std::cerr << "ERROR in LauAbsResonance : Total charge of daughters = " << totalCharge
127  << ". Resonance charge = " << resCharge_ << "." << std::endl;
128  gSystem->Exit( EXIT_FAILURE );
129  }
130 }
131 
132 // Constructor
133 LauAbsResonance::LauAbsResonance( const TString& resName,
134  const Int_t resPairAmpInt,
135  const LauDaughters* daughters,
136  const Int_t resSpin ) :
137  daughters_( daughters ),
138  resName_( resName ),
139  sanitisedName_( resName ),
140  resSpin_( resSpin ),
141  resPairAmpInt_( resPairAmpInt )
142 {
143  if ( daughters_ == 0 ) {
144  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided"
145  << std::endl;
146  gSystem->Exit( EXIT_FAILURE );
147  }
148 
149  nameParent_ = this->getNameParent();
150  nameDaug1_ = this->getNameDaug1();
151  nameDaug2_ = this->getNameDaug2();
152  nameBachelor_ = this->getNameBachelor();
153  massParent_ = this->getMassParent();
154  massDaug1_ = this->getMassDaug1();
155  massDaug2_ = this->getMassDaug2();
156  massBachelor_ = this->getMassBachelor();
157  chargeParent_ = this->getChargeParent();
158  chargeDaug1_ = this->getChargeDaug1();
159  chargeDaug2_ = this->getChargeDaug2();
161 
162  // Since we haven't been provided with a LauResonanceInfo object we can just
163  // set the change of the resonance to be the sum of the daughter charges
165 }
166 
167 // Destructor
169 {
170 }
171 
173 {
174  // Use LauKinematics interface for amplitude
175 
176  // For resonance made from tracks i, j, we need the momenta
177  // of tracks i and k in the i-j rest frame for spin helicity calculations
178  // in the Zemach tensor formalism.
179  // Also need the momentum of track k in the parent rest-frame for
180  // calculation of the Blatt-Weisskopf factors.
181  mass_ = 0.0;
182  cosHel_ = 0.0;
183  q_ = 0.0;
184  p_ = 0.0;
185  pstar_ = 0.0;
186  erm_ = 1.0;
187  covFactor_ = 1.0;
188 
189  if ( resPairAmpInt_ == 1 ) {
190 
191  mass_ = kinematics->getm23();
192  cosHel_ = kinematics->getc23();
193  q_ = kinematics->getp2_23();
194  p_ = kinematics->getp1_23();
195  pstar_ = kinematics->getp1_Parent();
196  erm_ = kinematics->getcov23();
197 
198  } else if ( resPairAmpInt_ == 2 ) {
199 
200  mass_ = kinematics->getm13();
201  cosHel_ = kinematics->getc13();
202  q_ = kinematics->getp1_13();
203  p_ = kinematics->getp2_13();
204  pstar_ = kinematics->getp2_Parent();
205  erm_ = kinematics->getcov13();
206 
207  } else if ( resPairAmpInt_ == 3 ) {
208 
209  mass_ = kinematics->getm12();
210  cosHel_ = kinematics->getc12();
211  q_ = kinematics->getp1_12();
212  p_ = kinematics->getp3_12();
213  pstar_ = kinematics->getp3_Parent();
214  erm_ = kinematics->getcov12();
215 
216  } else {
217  std::cerr << "ERROR in LauAbsResonance::amplitude : Nonsense setup of resPairAmp array."
218  << std::endl;
219  gSystem->Exit( EXIT_FAILURE );
220  }
221 
222  if ( this->flipHelicity() ) {
223  cosHel_ *= -1.0;
224  }
225 
226  if ( this->ignoreMomenta() ) {
227  q_ = 1.0;
228  p_ = 1.0;
229  pstar_ = 1.0;
230  erm_ = 1.0;
231  }
232 
233  // Calculate the spin factors
234  Double_t spinTerm( 1.0 );
235  Double_t pProd( 1.0 );
236 
237  if ( ! this->ignoreSpin() ) {
238  switch ( this->getSpinType() ) {
239 
240  case Zemach_P :
241  pProd = q_ * p_;
242  spinTerm = this->calcZemachSpinFactor( pProd );
243  break;
244 
245  case Zemach_Pstar :
246  pProd = q_ * pstar_;
247  spinTerm = this->calcZemachSpinFactor( pProd );
248  break;
249 
250  case Covariant :
251  pProd = q_ * pstar_;
252  spinTerm = this->calcCovSpinFactor( pProd );
253  break;
254 
255  case Legendre :
256  spinTerm = this->calcLegendrePoly();
257  break;
258  }
259  }
260 
261  // Calculate the full amplitude
262  LauComplex resAmplitude = this->resAmp( mass_, spinTerm );
263 
264  return resAmplitude;
265 }
266 
267 void LauAbsResonance::calcCovFactor( const Double_t erm )
268 {
269  if ( resSpin_ == 0 ) {
270  covFactor_ = 1.0;
271  } else if ( resSpin_ == 1 ) {
272  covFactor_ = erm;
273  } else if ( resSpin_ == 2 ) {
274  covFactor_ = erm * erm + 0.5;
275  } else if ( resSpin_ == 3 ) {
276  covFactor_ = erm * ( erm * erm + 1.5 );
277  } else if ( resSpin_ == 4 ) {
278  covFactor_ = ( 8. * erm * erm * erm * erm + 24. * erm * erm + 3. ) / 35.;
279  } else if ( resSpin_ > 4 ) {
280  std::cerr << "WARNING in LauAbsResonance::calcCovFactor : covariant spin factor cannot (yet) be fully calculated for spin >= 5"
281  << std::endl;
282  std::cerr << " : the function of sqrt(1 + (p/mParent)^2) part will be missing"
283  << std::endl;
284  covFactor_ = 1.0;
285  }
286 }
287 
288 Double_t LauAbsResonance::calcCovSpinFactor( const Double_t pProd )
289 {
290  if ( resSpin_ == 0 ) {
291  covFactor_ = 1.0;
292  return 1.0;
293  }
294 
295  // Covariant spin factor is (p* q)^L * f_L(erm) * P_L(cosHel)
296  Double_t spinFactor( pProd );
297  for ( Int_t i( 1 ); i < resSpin_; ++i ) {
298  spinFactor *= pProd;
299  }
300 
301  this->calcCovFactor( erm_ );
302 
303  spinFactor *= covFactor_;
304 
305  spinFactor *= this->calcLegendrePoly();
306 
307  return spinFactor;
308 }
309 
310 Double_t LauAbsResonance::calcZemachSpinFactor( const Double_t pProd ) const
311 {
312  // Calculate the spin factors
313  //
314  // These are calculated as follows
315  //
316  // -2^j * (q*p)^j * cj * Pj(cosHel)
317  //
318  // where Pj(coshHel) is the jth order Legendre polynomial and
319  //
320  // cj = j! / (2j-1)!!
321 
322  if ( resSpin_ == 0 ) {
323  return 1.0;
324  }
325 
326  Double_t spinFactor( pProd );
327  for ( Int_t i( 1 ); i < resSpin_; ++i ) {
328  spinFactor *= pProd;
329  }
330 
331  spinFactor *= this->calcLegendrePoly();
332 
333  return spinFactor;
334 }
335 
336 Double_t LauAbsResonance::calcLegendrePoly( const Double_t cosHel )
337 {
338  cosHel_ = cosHel;
339  return this->calcLegendrePoly();
340 }
341 
343 {
344  Double_t legPol = 1.0;
345 
346  if ( resSpin_ == 1 ) {
347  // Calculate vector resonance Legendre polynomial
348  legPol = -2.0 * cosHel_;
349  } else if ( resSpin_ == 2 ) {
350  // Calculate tensor resonance Legendre polynomial
351  legPol = 4.0 * ( 3.0 * cosHel_ * cosHel_ - 1.0 ) / 3.0;
352  } else if ( resSpin_ == 3 ) {
353  // Calculate spin 3 resonance Legendre polynomial
354  legPol = -8.0 * ( 5.0 * cosHel_ * cosHel_ * cosHel_ - 3.0 * cosHel_ ) / 5.0;
355  } else if ( resSpin_ == 4 ) {
356  // Calculate spin 4 resonance Legendre polynomial
357  legPol = 16.0 *
358  ( 35.0 * cosHel_ * cosHel_ * cosHel_ * cosHel_ - 30.0 * cosHel_ * cosHel_ + 3.0 ) /
359  35.0;
360  } else if ( resSpin_ == 5 ) {
361  // Calculate spin 5 resonance Legendre polynomial
362  legPol = -32.0 *
363  ( 63.0 * cosHel_ * cosHel_ * cosHel_ * cosHel_ * cosHel_ -
364  70.0 * cosHel_ * cosHel_ * cosHel_ + 15.0 * cosHel_ ) /
365  63.0;
366  } else if ( resSpin_ > 5 ) {
367  std::cerr << "WARNING in LauAbsResonance::calcLegendrePoly : Legendre polynomials not (yet) implemented for spin > 5"
368  << std::endl;
369  }
370 
371  return legPol;
372 }
373 
374 void LauAbsResonance::changeResonance( const Double_t newMass,
375  const Double_t newWidth,
376  const Int_t newSpin )
377 {
378  if ( newMass > 0.0 ) {
379  resMass_->valueAndRange( newMass, 0.0, 3.0 * newMass );
380  resMass_->initValue( newMass );
381  resMass_->genValue( newMass );
382  std::cout << "INFO in LauAbsResonance::changeResonance : Setting mass to "
383  << resMass_->value() << std::endl;
384  }
385  if ( newWidth > 0.0 ) {
386  resWidth_->valueAndRange( newWidth, 0.0, 3.0 * newWidth );
387  resWidth_->initValue( newWidth );
388  resWidth_->genValue( newWidth );
389  std::cout << "INFO in LauAbsResonance::changeResonance : Setting width to "
390  << resWidth_->value() << std::endl;
391  }
392  if ( newSpin > -1 ) {
393  resSpin_ = newSpin;
394  std::cout << "INFO in LauAbsResonance::changeResonance : Setting spin to " << resSpin_
395  << std::endl;
396  }
397 }
398 
399 void LauAbsResonance::changeBWBarrierRadii( const Double_t resRadius, const Double_t parRadius )
400 {
401  if ( resRadius >= 0.0 && resBWFactor_ != 0 ) {
402  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
403  resBWRadius->value( resRadius );
404  resBWRadius->initValue( resRadius );
405  resBWRadius->genValue( resRadius );
406  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting resonance factor radius to "
407  << resBWRadius->value() << std::endl;
408  }
409  if ( parRadius >= 0.0 && parBWFactor_ != 0 ) {
410  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
411  parBWRadius->value( parRadius );
412  parBWRadius->initValue( parRadius );
413  parBWRadius->genValue( parRadius );
414  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting parent factor radius to "
415  << parBWRadius->value() << std::endl;
416  }
417 }
418 
419 void LauAbsResonance::setResonanceParameter( const TString& name, const Double_t value )
420 {
421  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
422  std::cerr << "WARNING in LauAbsResonance::setResonanceParameter : Unable to set parameter \""
423  << name << "\" to value: " << value << "." << std::endl;
424 }
425 
427 {
428  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
429  std::cerr << "WARNING in LauAbsResonance::floatResonanceParameter : Unable to release parameter \""
430  << name << "\"." << std::endl;
431 }
432 
434 {
435  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
436  std::cerr << "WARNING in LauAbsResonance::getResonanceParameter : Unable to get parameter \""
437  << name << "\"." << std::endl;
438  return 0;
439 }
440 
442 {
443  if ( param == 0 ) {
444  return;
445  }
446 
447  if ( param->clone() ) {
448  resParameters_.push_back( param->parent() );
449  } else {
450  resParameters_.push_back( param );
451  }
452 }
453 
454 void LauAbsResonance::fixBarrierRadii( const Bool_t fixResRad, const Bool_t fixParRad )
455 {
456  if ( resBWFactor_ == 0 ) {
457  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : resonance barrier factor not present, cannot fix/float it"
458  << std::endl;
459  return;
460  }
461 
462  if ( parBWFactor_ == 0 ) {
463  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : parent barrier factor not present, cannot fix/float it"
464  << std::endl;
465  return;
466  }
467 
468  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
469  resBWRadius->fixed( fixResRad );
470 
471  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
472  parBWRadius->fixed( fixParRad );
473 }
474 
476 {
477  if ( resBWFactor_ == 0 ) {
478  std::cerr << "WARNING in LauAbsResonance::fixResRadius : resonance barrier factor not present"
479  << std::endl;
480  return kTRUE;
481  }
482 
484  return bwRadius->fixed();
485 }
486 
488 {
489  if ( parBWFactor_ == 0 ) {
490  std::cerr << "WARNING in LauAbsResonance::fixParRadius : parent barrier factor not present"
491  << std::endl;
492  return kTRUE;
493  }
494 
496  return bwRadius->fixed();
497 }
498 
500 {
501  if ( resBWFactor_ == 0 ) {
502  std::cerr << "WARNING in LauAbsResonance::getResRadius : resonance barrier factor not present"
503  << std::endl;
504  return -1.0;
505  }
506 
508  return bwRadius->unblindValue();
509 }
510 
512 {
513  if ( parBWFactor_ == 0 ) {
514  std::cerr << "WARNING in LauAbsResonance::getParRadius : parent barrier factor not present"
515  << std::endl;
516  return -1.0;
517  }
518 
520  return bwRadius->unblindValue();
521 }
522 
524 {
525  // Get the parent mass
526  Double_t mass( LauConstants::mB );
527 
528  if ( daughters_ ) {
529  mass = daughters_->getMassParent();
530  }
531 
532  return mass;
533 }
534 
536 {
537  // Get the daughter mass
538  Double_t mass( LauConstants::mPi );
539 
540  if ( daughters_ ) {
541  if ( resPairAmpInt_ == 1 ) {
542  mass = daughters_->getMassDaug2();
543  } else if ( resPairAmpInt_ == 2 ) {
544  mass = daughters_->getMassDaug1();
545  } else if ( resPairAmpInt_ == 3 ) {
546  mass = daughters_->getMassDaug1();
547  }
548  }
549 
550  return mass;
551 }
552 
554 {
555  // Get the daughter mass
556  Double_t mass( LauConstants::mPi );
557 
558  if ( daughters_ ) {
559  if ( resPairAmpInt_ == 1 ) {
560  mass = daughters_->getMassDaug3();
561  } else if ( resPairAmpInt_ == 2 ) {
562  mass = daughters_->getMassDaug3();
563  } else if ( resPairAmpInt_ == 3 ) {
564  mass = daughters_->getMassDaug2();
565  }
566  }
567 
568  return mass;
569 }
570 
572 {
573  // Get the bachelor mass
574  Double_t mass( LauConstants::mPi );
575 
576  if ( daughters_ ) {
577  if ( resPairAmpInt_ == 1 ) {
578  mass = daughters_->getMassDaug1();
579  } else if ( resPairAmpInt_ == 2 ) {
580  mass = daughters_->getMassDaug2();
581  } else if ( resPairAmpInt_ == 3 ) {
582  mass = daughters_->getMassDaug3();
583  }
584  }
585 
586  return mass;
587 }
588 
590 {
591  // Get the parent charge
592  Int_t charge( 0 );
593 
594  if ( daughters_ ) {
595  charge = daughters_->getChargeParent();
596  }
597 
598  return charge;
599 }
600 
602 {
603  // Get the daughter charge
604  Int_t charge( 0 );
605 
606  if ( daughters_ ) {
607  if ( resPairAmpInt_ == 1 ) {
608  charge = daughters_->getChargeDaug2();
609  } else if ( resPairAmpInt_ == 2 ) {
610  charge = daughters_->getChargeDaug1();
611  } else if ( resPairAmpInt_ == 3 ) {
612  charge = daughters_->getChargeDaug1();
613  }
614  }
615 
616  return charge;
617 }
618 
620 {
621  // Get the daughter charge
622  Int_t charge( 0 );
623 
624  if ( daughters_ ) {
625  if ( resPairAmpInt_ == 1 ) {
626  charge = daughters_->getChargeDaug3();
627  } else if ( resPairAmpInt_ == 2 ) {
628  charge = daughters_->getChargeDaug3();
629  } else if ( resPairAmpInt_ == 3 ) {
630  charge = daughters_->getChargeDaug2();
631  }
632  }
633 
634  return charge;
635 }
636 
638 {
639  // Get the bachelor charge
640  Int_t charge( 0 );
641 
642  if ( daughters_ ) {
643  if ( resPairAmpInt_ == 1 ) {
644  charge = daughters_->getChargeDaug1();
645  } else if ( resPairAmpInt_ == 2 ) {
646  charge = daughters_->getChargeDaug2();
647  } else if ( resPairAmpInt_ == 3 ) {
648  charge = daughters_->getChargeDaug3();
649  }
650  }
651 
652  return charge;
653 }
654 
656 {
657  // Get the parent name
658  TString name( "" );
659 
660  if ( daughters_ ) {
662  }
663 
664  return name;
665 }
666 
668 {
669  // Get the daughter name
670  TString name( "" );
671 
672  if ( daughters_ ) {
673  if ( resPairAmpInt_ == 1 ) {
675  } else if ( resPairAmpInt_ == 2 ) {
677  } else if ( resPairAmpInt_ == 3 ) {
679  }
680  }
681 
682  return name;
683 }
684 
686 {
687  // Get the daughter name
688  TString name( "" );
689 
690  if ( daughters_ ) {
691  if ( resPairAmpInt_ == 1 ) {
693  } else if ( resPairAmpInt_ == 2 ) {
695  } else if ( resPairAmpInt_ == 3 ) {
697  }
698  }
699 
700  return name;
701 }
702 
704 {
705  // Get the bachelor name
706  TString name( "" );
707 
708  if ( daughters_ ) {
709  if ( resPairAmpInt_ == 1 ) {
711  } else if ( resPairAmpInt_ == 2 ) {
713  } else if ( resPairAmpInt_ == 3 ) {
715  }
716  }
717 
718  return name;
719 }
Int_t resCharge_
Resonance charge.
TString getNameDaug3() const
Get name of the third daughter particle.
File containing declaration of LauResonanceInfo class.
Int_t getChargeDaug1() const
Get charge of the first daughter particle.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:49
Double_t unblindValue() const
The unblinded value of the parameter.
Double_t getc13() const
Get the cosine of the helicity angle theta13.
Double_t mass_
Invariant mass.
Int_t resSpin_
Resonance spin.
Int_t chargeBachelor_
Bachelor charge.
Double_t pstar_
Bachelor momentum in parent rest frame.
Double_t getp2_13() const
Get the momentum of the track 2 in 13 rest frame.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
TString getNameDaug2() const
Get name of the second daughter particle.
Double_t cosHel_
Helicity angle cosine.
Double_t value() const
The value of the parameter.
static bool isIncoherentModel(LauResonanceModel model)
Is the resonance model incoherent?
void calcCovFactor(const Double_t erm)
Calculate the spin-dependent covariant factor.
Double_t massDaug2_
Daughter 2 mass.
TString nameBachelor_
Bachelor name.
Double_t q_
Daughter momentum in resonance rest frame.
TString getNameParent() const
Get name of the parent particle.
Int_t chargeDaug2_
Daughter 2 charge.
Bool_t ignoreSpin() const
Get the ignore spin flag.
Double_t erm_
Covariant factor.
LauAbsResonance(LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor (for use by standard resonances)
LauBlattWeisskopfFactor * resBWFactor_
Blatt Weisskopf barrier for resonance decay.
LauParameter * resMass_
Resonance mass.
File containing declaration of LauParameter class.
File containing declaration of LauAbsResonance class.
Double_t massParent_
Parent mass.
void fixBarrierRadii(const Bool_t fixResRadius, const Bool_t fixParRadius)
Fix or release the Blatt-Weisskopf barrier radii.
Bool_t fixResRadius() const
Get the status of resonance barrier radius (fixed or released)
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)=0
Complex resonant amplitude.
Bool_t ignoreMomenta() const
Get the ignore momenta flag.
TString getNameDaug1() const
Get name of the first daughter particle.
Int_t getChargeParent() const
Get the Charge of the parent particle.
Double_t getp1_13() const
Get the momentum of the track 1 in 13 rest frame.
Class for defining a complex number.
Definition: LauComplex.hh:61
Double_t getc12() const
Get the cosine of the helicity angle theta12.
Double_t calcCovSpinFactor(const Double_t pProd)
Calculate the amplitude spin term using the covariant tensor formalism.
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 changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius)
Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed.
LauParameter * parent() const
The parent parameter.
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
LauParameter * resWidth_
Resonance width.
Int_t chargeParent_
Parent charge.
Int_t chargeDaug1_
Daughter 1 charge.
Double_t covFactor_
Covariant factor (full spin-dependent expression)
Bool_t clone() const
Check whether is a clone or not.
Double_t getm23() const
Get the m23 invariant mass.
Double_t getMassParent() const
Get the parent particle mass.
File containing declaration of LauDaughters class.
Double_t getp1_12() const
Get the momentum of the track 1 in 12 rest frame.
Double_t getp2_23() const
Get the momentum of the track 2 in 23 rest frame.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
Double_t getcov12() const
Get covariant factor in 12 axis.
Bool_t flipHelicity() const
Get the helicity flip flag.
TString getNameBachelor() const
Get the name of the daughter that does not originate form the resonance.
Double_t getm12() const
Get the m12 invariant mass.
Int_t getChargeDaug3() const
Get charge of the third daughter particle.
Double_t getMassParent() const
Get mass of the parent particle.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
Double_t getcov13() const
Get covariant factor in 13 axis.
Double_t getm13() const
Get the m13 invariant mass.
Class for defining the properties of a resonant particle.
virtual ~LauAbsResonance()
Destructor.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Calculate the complex amplitude.
Int_t getChargeDaug1() const
Get the charge of daughter 1.
Int_t resPairAmpInt_
DP axis identifier.
File containing LauConstants namespace.
const TString & name() const
The parameter name.
Double_t getp3_12() const
Get the momentum of the track 3 in 12 rest frame.
Double_t getMassDaug1() const
Get the mass of daughter 1.
TString nameDaug1_
Daughter 1 name.
TString getNameParent() const
Get the name of the parent particle.
Int_t getChargeDaug2() const
Get the charge of daughter 2.
const LauParameter * getRadiusParameter() const
Retrieve the radius parameter.
Int_t getChargeDaug2() const
Get charge of the second daughter particle.
Double_t getParRadius() const
Get the radius of the parent barrier factor.
Double_t getResRadius() const
Get the radius of the resonance barrier factor.
Int_t getChargeBachelor() const
Get the charge of the bachelor daughter.
std::vector< LauParameter * > resParameters_
All parameters of the resonance.
const LauDaughters * daughters_
Information on the particles.
Double_t getMassDaug2() const
Get the mass of daughter 2.
Bool_t fixParRadius() const
Get the status of parent barrier radius (fixed or released)
Double_t calcLegendrePoly() const
Calculate the Legendre polynomial for the spin factor.
Double_t calcZemachSpinFactor(const Double_t pProd) const
Calculate the amplitude spin term using the Zemach tensor formalism.
Double_t getp2_Parent() const
Get the momentum of the track 2 in parent rest frame.
Double_t massDaug1_
Daughter 1 mass.
Class for calculating 3-body kinematic quantities.
Double_t getp1_Parent() const
Get the momentum of the track 1 in parent rest frame.
Double_t genValue() const
The value generated for the parameter.
Double_t getp3_Parent() const
Get the momentum of the track 3 in parent rest frame.
Double_t massBachelor_
Bachelor mass.
TString nameParent_
Parent name.
const Double_t mPi
Mass of pi+- (GeV/c^2)
Definition: LauConstants.hh:55
Double_t getcov23() const
Get covariant factor in 23 axis.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getMassDaug2() const
Get mass of second daughter particle.
Double_t getMassDaug3() const
Get mass of third daughter particle.
const Double_t mB
Mass of charged B (GeV/c^2)
Definition: LauConstants.hh:49
TString getNameDaug1() const
Get the name of the first daughter of the resonance.
Int_t getChargeParent() const
Get charge of the parent particle.
void valueAndRange(Double_t newValue, Double_t newMinValue, Double_t newMaxValue)
Set the value and range for the parameter.
Double_t getMassBachelor() const
Get the mass of the bachelor daughter.
Double_t getMassDaug1() const
Get mass of first daughter particle.
TString getNameDaug2() const
Get the name of the second daughter of the resonance.
Double_t p_
Bachelor momentum in resonance rest frame.
LauSpinType getSpinType() const
Get the spin type.
TString nameDaug2_
Daughter 2 name.
Double_t initValue() const
The initial value of the parameter.
LauBlattWeisskopfFactor * parBWFactor_
Blatt Weisskopf barrier for parent decay.
LauResonanceModel
Define the allowed resonance types.
Double_t getp1_23() const
Get the momentum of the track 1 in 23 rest frame.
File containing declaration of LauKinematics class.