laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauAbsResonance.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2014.
3 // Distributed under the Boost Software License, Version 1.0.
4 // (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
5 
6 // Authors:
7 // Thomas Latham
8 // John Back
9 // Paul Harrison
10 
15 #include <iostream>
16 
17 #include "TSystem.h"
18 
19 #include "LauAbsResonance.hh"
20 #include "LauConstants.hh"
21 #include "LauDaughters.hh"
22 #include "LauKinematics.hh"
23 #include "LauParameter.hh"
24 #include "LauResonanceInfo.hh"
25 
27 
28 bool LauAbsResonance::isIncoherentModel(LauResonanceModel model) {
29  switch(model) {
30  case BW:
31  case RelBW:
32  case GS:
33  case Flatte:
34  case Sigma:
35  case Kappa:
36  case Dabba:
37  case LASS:
38  case LASS_BW:
39  case LASS_NR:
40  case EFKLLM:
41  case KMatrix:
42  case FlatNR:
43  case NRModel:
44  case BelleNR:
45  case PowerLawNR:
46  case BelleSymNR:
47  case BelleSymNRNoInter:
48  case TaylorNR:
49  case PolNR:
50  case MIPW_MagPhase:
51  case MIPW_RealImag:
52  case RhoOmegaMix_GS:
53  case RhoOmegaMix_RBW:
54  case RhoOmegaMix_GS_1:
55  case RhoOmegaMix_RBW_1:
56  break;
57  case GaussIncoh:
58  return true;
59  }
60  return false;
61 }
62 
63 // Constructor
64 LauAbsResonance::LauAbsResonance(LauResonanceInfo* resInfo, const Int_t resPairAmpInt, const LauDaughters* daughters) :
65  resInfo_(resInfo),
66  daughters_(daughters),
67  nameParent_(""), nameDaug1_(""), nameDaug2_(""), nameBachelor_(""),
68  chargeParent_(0), chargeDaug1_(0), chargeDaug2_(0), chargeBachelor_(0),
69  massParent_(0.0), massDaug1_(0.0), massDaug2_(0.0), massBachelor_(0.0),
70  resName_( (resInfo!=0) ? resInfo->getName() : "" ),
71  sanitisedName_( (resInfo!=0) ? resInfo->getSanitisedName() : "" ),
72  resMass_( (resInfo!=0) ? resInfo->getMass() : 0 ),
73  resWidth_( (resInfo!=0) ? resInfo->getWidth() : 0 ),
74  resSpin_( (resInfo!=0) ? resInfo->getSpin() : 0 ),
75  resCharge_( (resInfo!=0) ? resInfo->getCharge() : 0 ),
76  resPairAmpInt_(resPairAmpInt),
77  parBWFactor_(0),
78  resBWFactor_(0),
79  flipHelicity_(kFALSE),
80  ignoreMomenta_(kFALSE),
81  ignoreSpin_(kFALSE),
82  ignoreBarrierScaling_(kFALSE),
83  q_(0.0),
84  p_(0.0),
85  pstar_(0.0)
86 {
87  if ( resInfo == 0 ) {
88  std::cerr << "ERROR in LauAbsResonance constructor : null LauResonanceInfo object provided" << std::endl;
89  gSystem->Exit(EXIT_FAILURE);
90  }
91 
92  if ( daughters_ == 0 ) {
93  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl;
94  gSystem->Exit(EXIT_FAILURE);
95  }
96 
97  nameParent_ = this->getNameParent();
98  nameDaug1_ = this->getNameDaug1();
99  nameDaug2_ = this->getNameDaug2();
100  nameBachelor_ = this->getNameBachelor();
101  massParent_ = this->getMassParent();
102  massDaug1_ = this->getMassDaug1();
103  massDaug2_ = this->getMassDaug2();
104  massBachelor_ = this->getMassBachelor();
105  chargeParent_ = this->getChargeParent();
106  chargeDaug1_ = this->getChargeDaug1();
107  chargeDaug2_ = this->getChargeDaug2();
109 
110  // check that the total charge adds up to that of the resonance
111  Int_t totalCharge = chargeDaug1_ + chargeDaug2_;
112  if ( (totalCharge != resCharge_) && (resPairAmpInt_ != 0) ) {
113  std::cerr << "ERROR in LauAbsResonance : Total charge of daughters = " << totalCharge << ". Resonance charge = " << resCharge_ << "." << std::endl;
114  gSystem->Exit(EXIT_FAILURE);
115  }
116 }
117 
118 // Constructor
119 LauAbsResonance::LauAbsResonance(const TString& resName, const Int_t resPairAmpInt, const LauDaughters* daughters) :
120  resInfo_(0),
121  daughters_(daughters),
122  nameParent_(""), nameDaug1_(""), nameDaug2_(""), nameBachelor_(""),
123  chargeParent_(0), chargeDaug1_(0), chargeDaug2_(0), chargeBachelor_(0),
124  massParent_(0.0), massDaug1_(0.0), massDaug2_(0.0), massBachelor_(0.0),
125  resName_(resName),
126  sanitisedName_(resName),
127  resMass_(0),
128  resWidth_(0),
129  resSpin_(0),
130  resCharge_(0),
131  resPairAmpInt_(resPairAmpInt),
132  parBWFactor_(0),
133  resBWFactor_(0),
134  flipHelicity_(kFALSE),
135  ignoreMomenta_(kFALSE),
136  ignoreSpin_(kFALSE),
137  ignoreBarrierScaling_(kFALSE),
138  q_(0.0),
139  p_(0.0),
140  pstar_(0.0)
141 {
142  if ( daughters_ == 0 ) {
143  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl;
144  gSystem->Exit(EXIT_FAILURE);
145  }
146 
147  nameParent_ = this->getNameParent();
148  nameDaug1_ = this->getNameDaug1();
149  nameDaug2_ = this->getNameDaug2();
150  nameBachelor_ = this->getNameBachelor();
151  massParent_ = this->getMassParent();
152  massDaug1_ = this->getMassDaug1();
153  massDaug2_ = this->getMassDaug2();
154  massBachelor_ = this->getMassBachelor();
155  chargeParent_ = this->getChargeParent();
156  chargeDaug1_ = this->getChargeDaug1();
157  chargeDaug2_ = this->getChargeDaug2();
159 
160  // check that the total charge adds up to that of the resonance
161  Int_t totalCharge = chargeDaug1_ + chargeDaug2_;
162  if ( (totalCharge != resCharge_) && (resPairAmpInt_ != 0) ) {
163  std::cerr << "ERROR in LauAbsResonance : Total charge of daughters = " << totalCharge << ". Resonance charge = " << resCharge_ << "." << std::endl;
164  gSystem->Exit(EXIT_FAILURE);
165  }
166 }
167 
168 // Destructor
170 {
171 }
172 
174 {
175  // Use LauKinematics interface for amplitude
176  Double_t mass(0.0), cosHel(0.0);
177  // For resonance made from tracks i, j, we need the momenta
178  // of tracks i and k in the i-j rest frame for spin helicity calculations
179  // in the Zemach tensor formalism.
180  // Also need the momentum of track k in the parent rest-frame for
181  // calculation of the Blatt-Weisskopf factors.
182  q_ = 0.0; p_ = 0.0; pstar_ = 0.0;
183 
184  if (resPairAmpInt_ == 1) {
185 
186  mass = kinematics->getm23();
187  cosHel = kinematics->getc23();
188  q_ = kinematics->getp2_23();
189  p_ = kinematics->getp1_23();
190  pstar_ = kinematics->getp1_Parent();
191 
192  } else if (resPairAmpInt_ == 2) {
193 
194  mass = kinematics->getm13();
195  cosHel = kinematics->getc13();
196  q_ = kinematics->getp1_13();
197  p_ = kinematics->getp2_13();
198  pstar_ = kinematics->getp2_Parent();
199 
200  } else if (resPairAmpInt_ == 3) {
201 
202  mass = kinematics->getm12();
203  cosHel = kinematics->getc12();
204  q_ = kinematics->getp1_12();
205  p_ = kinematics->getp3_12();
206  pstar_ = kinematics->getp3_Parent();
207 
208  } else {
209  std::cerr << "ERROR in LauAbsResonance::amplitude : Nonsense setup of resPairAmp array." << std::endl;
210  gSystem->Exit(EXIT_FAILURE);
211  }
212 
213  if (this->flipHelicity()) {
214  cosHel *= -1.0;
215  }
216 
217  if (this->ignoreMomenta()) {
218  q_ = 1.0;
219  p_ = 1.0;
220  }
221 
222  // Calculate the spin factors
223  Double_t spinTerm(1.0);
224  if (!this->ignoreSpin()) {
225  Double_t pProd = q_*p_;
226  spinTerm = this->calcSpinTerm( cosHel, pProd );
227  }
228 
229  // Calculate the full amplitude
230  LauComplex resAmplitude = this->resAmp(mass, spinTerm);
231 
232  return resAmplitude;
233 }
234 
235 Double_t LauAbsResonance::calcSpinTerm( const Double_t cosHel, const Double_t pProd ) const
236 {
237  // Calculate the spin factors
238  //
239  // These are calculated as follows
240  //
241  // -2^j * (q*p)^j * cj * Pj(cosHel)
242  //
243  // where Pj(coshHel) is the jth order Legendre polynomial and
244  //
245  // cj = j! / (2j-1)!!
246 
247  Double_t spinTerm = 1.0;
248 
249  if (resSpin_ == 1) {
250  // Calculate vector resonance Zemach helicity factor
251  spinTerm = -2.0*pProd*cosHel;
252  } else if (resSpin_ == 2) {
253  // Calculate tensor resonance Zemach helicity factor
254  spinTerm = 4.0*(pProd*pProd)*(3.0*cosHel*cosHel - 1.0)/3.0;
255  } else if (resSpin_ == 3) {
256  // Calculate spin 3 resonance Zemach helicity factor
257  spinTerm = -8.0*(pProd*pProd*pProd)*(5.0*cosHel*cosHel*cosHel - 3.0*cosHel)/5.0;
258  } else if (resSpin_ == 4) {
259  // Calculate spin 4 resonance Zemach helicity factor
260  spinTerm = 16.0*(pProd*pProd*pProd*pProd)*(35.0*cosHel*cosHel*cosHel*cosHel - 30.0*cosHel*cosHel + 3.0)/35.0;
261  } else if (resSpin_ == 5) {
262  // Calculate spin 5 resonance Zemach helicity factor
263  spinTerm = -32.0*(pProd*pProd*pProd*pProd*pProd)*(63.0*cosHel*cosHel*cosHel*cosHel*cosHel - 70.0*cosHel*cosHel*cosHel + 15.0*cosHel)/63.0;
264  }
265 
266  return spinTerm;
267 }
268 
269 void LauAbsResonance::changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin)
270 {
271  if (newMass > 0.0) {
272  resMass_->valueAndRange(newMass,0.0,3.0*newMass);
273  resMass_->initValue(newMass);
274  resMass_->genValue(newMass);
275  std::cout << "INFO in LauAbsResonance::changeResonance : Setting mass to " << resMass_->value() << std::endl;
276  }
277  if (newWidth > 0.0) {
278  resWidth_->valueAndRange(newWidth,0.0,3.0*newWidth);
279  resWidth_->initValue(newWidth);
280  resWidth_->genValue(newWidth);
281  std::cout << "INFO in LauAbsResonance::changeResonance : Setting width to " << resWidth_->value() << std::endl;
282  }
283  if (newSpin > -1) {
284  resSpin_ = newSpin;
285  std::cout << "INFO in LauAbsResonance::changeResonance : Setting spin to " << resSpin_ << std::endl;
286  }
287 }
288 
289 void LauAbsResonance::changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius)
290 {
291  if ( resRadius >= 0.0 && resBWFactor_ != 0 ) {
292  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
293  resBWRadius->value(resRadius);
294  resBWRadius->initValue(resRadius);
295  resBWRadius->genValue(resRadius);
296  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting resonance factor radius to " << resBWRadius->value() << std::endl;
297  }
298  if ( parRadius >= 0.0 && parBWFactor_ != 0 ) {
299  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
300  parBWRadius->value(parRadius);
301  parBWRadius->initValue(parRadius);
302  parBWRadius->genValue(parRadius);
303  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting parent factor radius to " << parBWRadius->value() << std::endl;
304  }
305 }
306 
307 void LauAbsResonance::setResonanceParameter(const TString& name, const Double_t value)
308 {
309  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
310  std::cerr << "WARNING in LauAbsResonance::setResonanceParameter : Unable to set parameter \"" << name << "\" to value: " << value << "." << std::endl;
311 }
312 
314 {
315  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
316  std::cerr << "WARNING in LauAbsResonance::floatResonanceParameter : Unable to release parameter \"" << name << "\"." << std::endl;
317 }
318 
320 {
321  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
322  std::cerr << "WARNING in LauAbsResonance::getResonanceParameter : Unable to get parameter \"" << name << "\"." << std::endl;
323  return 0;
324 }
325 
327 {
328  if ( param == 0 ) {
329  return;
330  }
331 
332  if ( param->clone() ) {
333  resParameters_.push_back( param->parent() );
334  } else {
335  resParameters_.push_back( param );
336  }
337 }
338 
339 void LauAbsResonance::fixBarrierRadii(const Bool_t fixResRad, const Bool_t fixParRad)
340 {
341  if ( resBWFactor_ == 0 ) {
342  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : resonance barrier factor not present, cannot fix/float it" << std::endl;
343  return;
344  }
345 
346  if ( parBWFactor_ == 0 ) {
347  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : parent barrier factor not present, cannot fix/float it" << std::endl;
348  return;
349  }
350 
351  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
352  resBWRadius->fixed(fixResRad);
353 
354  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
355  parBWRadius->fixed(fixParRad);
356 }
357 
359 {
360  if ( resBWFactor_ == 0 ) {
361  std::cerr << "WARNING in LauAbsResonance::fixResRadius : resonance barrier factor not present" << std::endl;
362  return kTRUE;
363  }
364 
366  return bwRadius->fixed();
367 }
368 
370 {
371  if ( parBWFactor_ == 0 ) {
372  std::cerr << "WARNING in LauAbsResonance::fixParRadius : parent barrier factor not present" << std::endl;
373  return kTRUE;
374  }
375 
377  return bwRadius->fixed();
378 }
379 
381 {
382  if ( resBWFactor_ == 0 ) {
383  std::cerr << "WARNING in LauAbsResonance::getResRadius : resonance barrier factor not present" << std::endl;
384  return -1.0;
385  }
386 
388  return bwRadius->unblindValue();
389 }
390 
392 {
393  if ( parBWFactor_ == 0 ) {
394  std::cerr << "WARNING in LauAbsResonance::getParRadius : parent barrier factor not present" << std::endl;
395  return -1.0;
396  }
397 
399  return bwRadius->unblindValue();
400 }
401 
403 {
404  // Get the parent mass
405  Double_t mass(LauConstants::mB);
406 
407  if (daughters_) {
408  mass = daughters_->getMassParent();
409  }
410 
411  return mass;
412 }
413 
415 {
416  // Get the daughter mass
417  Double_t mass(LauConstants::mPi);
418 
419  if (daughters_) {
420  if (resPairAmpInt_ == 1) {
421  mass = daughters_->getMassDaug2();
422  } else if (resPairAmpInt_ == 2) {
423  mass = daughters_->getMassDaug1();
424  } else if (resPairAmpInt_ == 3) {
425  mass = daughters_->getMassDaug1();
426  }
427  }
428 
429  return mass;
430 }
431 
433 {
434  // Get the daughter mass
435  Double_t mass(LauConstants::mPi);
436 
437  if (daughters_) {
438  if (resPairAmpInt_ == 1) {
439  mass = daughters_->getMassDaug3();
440  } else if (resPairAmpInt_ == 2) {
441  mass = daughters_->getMassDaug3();
442  } else if (resPairAmpInt_ == 3) {
443  mass = daughters_->getMassDaug2();
444  }
445  }
446 
447  return mass;
448 }
449 
451 {
452  // Get the bachelor mass
453  Double_t mass(LauConstants::mPi);
454 
455  if (daughters_) {
456  if (resPairAmpInt_ == 1) {
457  mass = daughters_->getMassDaug1();
458  } else if (resPairAmpInt_ == 2) {
459  mass = daughters_->getMassDaug2();
460  } else if (resPairAmpInt_ == 3) {
461  mass = daughters_->getMassDaug3();
462  }
463  }
464 
465  return mass;
466 }
467 
469 {
470  // Get the parent charge
471  Int_t charge(0);
472 
473  if (daughters_) {
474  charge = daughters_->getChargeParent();
475  }
476 
477  return charge;
478 }
479 
481 {
482  // Get the daughter charge
483  Int_t charge(0);
484 
485  if (daughters_) {
486  if (resPairAmpInt_ == 1) {
487  charge = daughters_->getChargeDaug2();
488  } else if (resPairAmpInt_ == 2) {
489  charge = daughters_->getChargeDaug1();
490  } else if (resPairAmpInt_ == 3) {
491  charge = daughters_->getChargeDaug1();
492  }
493  }
494 
495  return charge;
496 }
497 
499 {
500  // Get the daughter charge
501  Int_t charge(0);
502 
503  if (daughters_) {
504  if (resPairAmpInt_ == 1) {
505  charge = daughters_->getChargeDaug3();
506  } else if (resPairAmpInt_ == 2) {
507  charge = daughters_->getChargeDaug3();
508  } else if (resPairAmpInt_ == 3) {
509  charge = daughters_->getChargeDaug2();
510  }
511  }
512 
513  return charge;
514 }
515 
517 {
518  // Get the bachelor charge
519  Int_t charge(0);
520 
521  if (daughters_) {
522  if (resPairAmpInt_ == 1) {
523  charge = daughters_->getChargeDaug1();
524  } else if (resPairAmpInt_ == 2) {
525  charge = daughters_->getChargeDaug2();
526  } else if (resPairAmpInt_ == 3) {
527  charge = daughters_->getChargeDaug3();
528  }
529  }
530 
531  return charge;
532 }
533 
535 {
536  // Get the parent name
537  TString name("");
538 
539  if (daughters_) {
540  name = daughters_->getNameParent();
541  }
542 
543  return name;
544 }
545 
547 {
548  // Get the daughter name
549  TString name("");
550 
551  if (daughters_) {
552  if (resPairAmpInt_ == 1) {
553  name = daughters_->getNameDaug2();
554  } else if (resPairAmpInt_ == 2) {
555  name = daughters_->getNameDaug1();
556  } else if (resPairAmpInt_ == 3) {
557  name = daughters_->getNameDaug1();
558  }
559  }
560 
561  return name;
562 }
563 
565 {
566  // Get the daughter name
567  TString name("");
568 
569  if (daughters_) {
570  if (resPairAmpInt_ == 1) {
571  name = daughters_->getNameDaug3();
572  } else if (resPairAmpInt_ == 2) {
573  name = daughters_->getNameDaug3();
574  } else if (resPairAmpInt_ == 3) {
575  name = daughters_->getNameDaug2();
576  }
577  }
578 
579  return name;
580 }
581 
583 {
584  // Get the bachelor name
585  TString name("");
586 
587  if (daughters_) {
588  if (resPairAmpInt_ == 1) {
589  name = daughters_->getNameDaug1();
590  } else if (resPairAmpInt_ == 2) {
591  name = daughters_->getNameDaug2();
592  } else if (resPairAmpInt_ == 3) {
593  name = daughters_->getNameDaug3();
594  }
595  }
596 
597  return name;
598 }
599 
Double_t getMassParent() const
Get mass of the parent particle.
TString getNameDaug1() const
Get name of the first daughter particle.
Int_t getChargeDaug3() const
Get charge of the third daughter particle.
Int_t getChargeBachelor() const
Get the charge of the bachelor daughter.
Double_t getc23() const
Get the cosine of the helicity angle theta23.
Bool_t fixed() const
Check whether the parameter is fixed or floated.
Double_t massDaug2_
Daughter 2 mass.
Double_t getMassBachelor() const
Get the mass of the bachelor daughter.
LauAbsResonance(LauResonanceInfo *resInfo, const Int_t resPairAmpInt, const LauDaughters *daughters)
Constructor (for use by standard resonances)
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.
Double_t getp2_Parent() const
Get the momentum of the track 2 in parent rest frame.
File containing declaration of LauResonanceInfo class.
ClassImp(LauAbsCoeffSet)
Double_t getMassParent() const
Get the parent particle mass.
Double_t getm23() const
Get the m23 invariant mass.
Double_t massDaug1_
Daughter 1 mass.
Class for defining the properties of a resonant particle.
const TString & name() const
The parameter name.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:33
Double_t getc13() const
Get the cosine of the helicity angle theta13.
TString nameDaug2_
Daughter 2 name.
Double_t getp2_13() const
Get the momentum of the track 2 in 13 rest frame.
File containing declaration of LauDaughters class.
virtual void floatResonanceParameter(const TString &name)
Allow the various parameters to float in the fit.
Int_t getChargeDaug1() const
Get charge of the first daughter particle.
const LauParameter * getRadiusParameter() const
Retrieve the radius parameter.
TString getNameDaug2() const
Get the name of the second daughter of the resonance.
LauParameter * resWidth_
Resonance width.
LauParameter * resMass_
Resonance mass.
LauBlattWeisskopfFactor * resBWFactor_
Blatt Weisskopf barrier for resonance decay.
Double_t getp1_12() const
Get the momentum of the track 1 in 12 rest frame.
Double_t massParent_
Parent mass.
Double_t getMassDaug1() const
Get the mass of daughter 1.
TString getNameDaug2() const
Get name of the second daughter particle.
Double_t getMassDaug2() const
Get mass of second daughter particle.
void valueAndRange(Double_t newValue, Double_t newMinValue, Double_t newMaxValue)
Set the value and range for the parameter.
Bool_t fixResRadius() const
Get the status of resonance barrier radius (fixed or released)
Bool_t flipHelicity() const
Get the helicity flip flag.
Double_t getParRadius() const
Get the radius of the parent barrier factor.
File containing declaration of LauKinematics class.
Double_t getMassDaug2() const
Get the mass of daughter 2.
virtual LauParameter * getResonanceParameter(const TString &name)
Access the given resonance parameter.
Double_t getp1_Parent() const
Get the momentum of the track 1 in parent rest frame.
const Double_t mPi
Mass of pi+- (GeV/c^2)
Definition: LauConstants.hh:40
void addFloatingParameter(LauParameter *param)
Add parameter to the list of floating parameters.
Bool_t clone() const
Check whether is a clone or not.
const Double_t mB
Mass of charged B (GeV/c^2)
Definition: LauConstants.hh:34
Int_t chargeBachelor_
Bachelor charge.
Int_t getChargeDaug1() const
Get the charge of daughter 1.
const LauDaughters * daughters_
Information on the particles.
Double_t getp1_23() const
Get the momentum of the track 1 in 23 rest frame.
File containing declaration of LauParameter class.
std::vector< LauParameter * > resParameters_
All parameters of the resonance.
Double_t getc12() const
Get the cosine of the helicity angle theta12.
Double_t getm13() const
Get the m13 invariant mass.
Int_t getChargeDaug2() const
Get charge of the second daughter particle.
virtual void setResonanceParameter(const TString &name, const Double_t value)
Set value of the various parameters.
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.
virtual ~LauAbsResonance()
Destructor.
LauBlattWeisskopfFactor * parBWFactor_
Blatt Weisskopf barrier for parent decay.
virtual LauComplex resAmp(Double_t mass, Double_t spinTerm)=0
Complex resonant amplitude.
TString getNameParent() const
Get the name of the parent particle.
TString nameBachelor_
Bachelor name.
Int_t chargeDaug2_
Daughter 2 charge.
Class for defining the fit parameter objects.
Definition: LauParameter.hh:35
Bool_t ignoreMomenta() const
Get the ignore momenta flag.
Int_t resCharge_
Resonance charge.
TString nameParent_
Parent name.
LauParameter * parent() const
The parent parameter.
Int_t chargeParent_
Parent charge.
TString getNameBachelor() const
Get the name of the daughter that does not originate form the resonance.
Int_t getChargeDaug2() const
Get the charge of daughter 2.
LauResonanceModel
Define the allowed resonance types.
Double_t getResRadius() const
Get the radius of the resonance barrier factor.
Double_t calcSpinTerm(const Double_t cosHel, const Double_t pProd) const
Calculate the amplitude spin term.
Double_t getp3_12() const
Get the momentum of the track 3 in 12 rest frame.
Bool_t fixParRadius() const
Get the status of parent barrier radius (fixed or released)
void fixBarrierRadii(const Bool_t fixResRadius, const Bool_t fixParRadius)
Fix or release the Blatt-Weisskopf barrier radii.
Abstract class for defining type for resonance amplitude models (Breit-Wigner, Flatte etc...
TString getNameParent() const
Get name of the parent particle.
void changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius)
Allow the Blatt-Weisskopf radius for the resonance and parent factors to be changed.
Int_t getChargeParent() const
Get charge of the parent particle.
Double_t initValue() const
The initial value of the parameter.
Double_t getMassDaug1() const
Get mass of first daughter particle.
Int_t resPairAmpInt_
DP axis identifier.
Double_t getm12() const
Get the m12 invariant mass.
File containing LauConstants namespace.
File containing declaration of LauAbsResonance class.
Double_t getp2_23() const
Get the momentum of the track 2 in 23 rest frame.
Double_t unblindValue() const
The unblinded value of the parameter.
Int_t chargeDaug1_
Daughter 1 charge.
Class for defining a complex number.
Definition: LauComplex.hh:47
TString getNameDaug1() const
Get the name of the first daughter of the resonance.
Class for calculating 3-body kinematic quantities.
Int_t resSpin_
Resonance spin.
Double_t value() const
The value of the parameter.
Bool_t ignoreSpin() const
Get the ignore spin flag.
TString nameDaug1_
Daughter 1 name.
virtual LauComplex amplitude(const LauKinematics *kinematics)
Calculate the complex amplitude.
Double_t getMassDaug3() const
Get mass of third daughter particle.
TString getNameDaug3() const
Get name of the third daughter particle.
Double_t genValue() const
The value generated for the parameter.
Double_t p_
Bachelor momentum in resonance rest frame.
Double_t getp3_Parent() const
Get the momentum of the track 3 in parent rest frame.
Double_t pstar_
Bachelor momentum in parent rest frame.
Double_t q_
Daughter momentum in resonance rest frame.