laura is hosted by Hepforge, IPPP Durham
Laura++  v3r3
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  spinType_(Zemach_P),
80  flipHelicity_(kFALSE),
81  ignoreMomenta_(kFALSE),
82  ignoreSpin_(kFALSE),
83  ignoreBarrierScaling_(kFALSE),
84  mass_(0.0),
85  cosHel_(0.0),
86  q_(0.0),
87  p_(0.0),
88  pstar_(0.0),
89  erm_(1.0),
90  covFactor_(1.0)
91 {
92  if ( resInfo == 0 ) {
93  std::cerr << "ERROR in LauAbsResonance constructor : null LauResonanceInfo object provided" << std::endl;
94  gSystem->Exit(EXIT_FAILURE);
95  }
96 
97  if ( daughters_ == 0 ) {
98  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl;
99  gSystem->Exit(EXIT_FAILURE);
100  }
101 
102  nameParent_ = this->getNameParent();
103  nameDaug1_ = this->getNameDaug1();
104  nameDaug2_ = this->getNameDaug2();
105  nameBachelor_ = this->getNameBachelor();
106  massParent_ = this->getMassParent();
107  massDaug1_ = this->getMassDaug1();
108  massDaug2_ = this->getMassDaug2();
109  massBachelor_ = this->getMassBachelor();
110  chargeParent_ = this->getChargeParent();
111  chargeDaug1_ = this->getChargeDaug1();
112  chargeDaug2_ = this->getChargeDaug2();
114 
115  // check that the total charge adds up to that of the resonance
116  Int_t totalCharge = chargeDaug1_ + chargeDaug2_;
117  if ( (totalCharge != resCharge_) && (resPairAmpInt_ != 0) ) {
118  std::cerr << "ERROR in LauAbsResonance : Total charge of daughters = " << totalCharge << ". Resonance charge = " << resCharge_ << "." << std::endl;
119  gSystem->Exit(EXIT_FAILURE);
120  }
121 }
122 
123 // Constructor
124 LauAbsResonance::LauAbsResonance(const TString& resName, const Int_t resPairAmpInt, const LauDaughters* daughters) :
125  resInfo_(0),
126  daughters_(daughters),
127  nameParent_(""), nameDaug1_(""), nameDaug2_(""), nameBachelor_(""),
128  chargeParent_(0), chargeDaug1_(0), chargeDaug2_(0), chargeBachelor_(0),
129  massParent_(0.0), massDaug1_(0.0), massDaug2_(0.0), massBachelor_(0.0),
130  resName_(resName),
131  sanitisedName_(resName),
132  resMass_(0),
133  resWidth_(0),
134  resSpin_(0),
135  resCharge_(0),
136  resPairAmpInt_(resPairAmpInt),
137  parBWFactor_(0),
138  resBWFactor_(0),
139  spinType_(Zemach_P),
140  flipHelicity_(kFALSE),
141  ignoreMomenta_(kFALSE),
142  ignoreSpin_(kFALSE),
143  ignoreBarrierScaling_(kFALSE),
144  mass_(0.0),
145  cosHel_(0.0),
146  q_(0.0),
147  p_(0.0),
148  pstar_(0.0),
149  erm_(1.0),
150  covFactor_(1.0)
151 {
152  if ( daughters_ == 0 ) {
153  std::cerr << "ERROR in LauAbsResonance constructor : null LauDaughters object provided" << std::endl;
154  gSystem->Exit(EXIT_FAILURE);
155  }
156 
157  nameParent_ = this->getNameParent();
158  nameDaug1_ = this->getNameDaug1();
159  nameDaug2_ = this->getNameDaug2();
160  nameBachelor_ = this->getNameBachelor();
161  massParent_ = this->getMassParent();
162  massDaug1_ = this->getMassDaug1();
163  massDaug2_ = this->getMassDaug2();
164  massBachelor_ = this->getMassBachelor();
165  chargeParent_ = this->getChargeParent();
166  chargeDaug1_ = this->getChargeDaug1();
167  chargeDaug2_ = this->getChargeDaug2();
169 
170  // Since we haven't been provided with a LauResonanceInfo object we can just
171  // set the change of the resonance to be the sum of the daughter charges
173 }
174 
175 // Destructor
177 {
178 }
179 
181 {
182  // Use LauKinematics interface for amplitude
183 
184  // For resonance made from tracks i, j, we need the momenta
185  // of tracks i and k in the i-j rest frame for spin helicity calculations
186  // in the Zemach tensor formalism.
187  // Also need the momentum of track k in the parent rest-frame for
188  // calculation of the Blatt-Weisskopf factors.
189  mass_ = 0.0; cosHel_ = 0.0;
190  q_ = 0.0; p_ = 0.0; pstar_ = 0.0;
191  erm_ = 1.0; covFactor_ = 1.0;
192 
193  if (resPairAmpInt_ == 1) {
194 
195  mass_ = kinematics->getm23();
196  cosHel_ = kinematics->getc23();
197  q_ = kinematics->getp2_23();
198  p_ = kinematics->getp1_23();
199  pstar_ = kinematics->getp1_Parent();
200  erm_ = kinematics->getcov23();
201 
202  } else if (resPairAmpInt_ == 2) {
203 
204  mass_ = kinematics->getm13();
205  cosHel_ = kinematics->getc13();
206  q_ = kinematics->getp1_13();
207  p_ = kinematics->getp2_13();
208  pstar_ = kinematics->getp2_Parent();
209  erm_ = kinematics->getcov13();
210 
211  } else if (resPairAmpInt_ == 3) {
212 
213  mass_ = kinematics->getm12();
214  cosHel_ = kinematics->getc12();
215  q_ = kinematics->getp1_12();
216  p_ = kinematics->getp3_12();
217  pstar_ = kinematics->getp3_Parent();
218  erm_ = kinematics->getcov12();
219 
220  } else {
221  std::cerr << "ERROR in LauAbsResonance::amplitude : Nonsense setup of resPairAmp array." << std::endl;
222  gSystem->Exit(EXIT_FAILURE);
223  }
224 
225  if (this->flipHelicity()) {
226  cosHel_ *= -1.0;
227  }
228 
229  if (this->ignoreMomenta()) {
230  q_ = 1.0;
231  p_ = 1.0;
232  pstar_ = 1.0;
233  erm_ = 1.0;
234  }
235 
236  // Calculate the spin factors
237  Double_t spinTerm(1.0);
238  Double_t pProd(1.0);
239 
240  if (!this->ignoreSpin()) {
241  switch ( this->getSpinType() ) {
242 
243  case Zemach_P:
244  pProd = q_*p_;
245  spinTerm = this->calcZemachSpinFactor( pProd );
246  break;
247 
248  case Zemach_Pstar:
249  pProd = q_*pstar_;
250  spinTerm = this->calcZemachSpinFactor( pProd );
251  break;
252 
253  case Covariant:
254  pProd = q_*pstar_;
255  spinTerm = this->calcCovSpinFactor( pProd );
256  break;
257 
258  case Legendre:
259  spinTerm = this->calcLegendrePoly();
260  break;
261  }
262  }
263 
264  // Calculate the full amplitude
265  LauComplex resAmplitude = this->resAmp(mass_, spinTerm);
266 
267  return resAmplitude;
268 }
269 
270 void LauAbsResonance::calcCovFactor( const Double_t erm )
271 {
272  if (resSpin_ == 0) {
273  covFactor_ = 1.0;
274  } else if (resSpin_ == 1) {
275  covFactor_ = erm;
276  } else if (resSpin_ == 2) {
277  covFactor_ = erm*erm + 0.5;
278  } else if (resSpin_ == 3) {
279  covFactor_ = erm*(erm*erm + 1.5);
280  } else if (resSpin_ == 4) {
281  covFactor_ = (8.*erm*erm*erm*erm + 24.*erm*erm + 3.)/35.;
282  } else if (resSpin_ > 4) {
283  std::cerr << "WARNING in LauAbsResonance::calcCovFactor : covariant spin factor cannot (yet) be fully calculated for spin >= 5" << std::endl;
284  std::cerr << " : the function of sqrt(1 + (p/mParent)^2) part will be missing" << std::endl;
285  covFactor_ = 1.0;
286  }
287 }
288 
289 Double_t LauAbsResonance::calcCovSpinFactor( const Double_t pProd )
290 {
291  if (resSpin_ == 0) {
292  covFactor_ = 1.0;
293  return 1.0;
294  }
295 
296  // Covariant spin factor is (p* q)^L * f_L(erm) * P_L(cosHel)
297  Double_t spinFactor(pProd);
298  for ( Int_t i(1); i < resSpin_; ++i ) {
299  spinFactor *= pProd;
300  }
301 
302  this->calcCovFactor( erm_ );
303 
304  spinFactor *= covFactor_;
305 
306  spinFactor *= this->calcLegendrePoly();
307 
308  return spinFactor;
309 }
310 
311 Double_t LauAbsResonance::calcZemachSpinFactor( const Double_t pProd ) const
312 {
313  // Calculate the spin factors
314  //
315  // These are calculated as follows
316  //
317  // -2^j * (q*p)^j * cj * Pj(cosHel)
318  //
319  // where Pj(coshHel) is the jth order Legendre polynomial and
320  //
321  // cj = j! / (2j-1)!!
322 
323  if (resSpin_ == 0) {
324  return 1.0;
325  }
326 
327  Double_t spinFactor(pProd);
328  for ( Int_t i(1); i < resSpin_; ++i ) {
329  spinFactor *= pProd;
330  }
331 
332  spinFactor *= this->calcLegendrePoly();
333 
334  return spinFactor;
335 }
336 
337 Double_t LauAbsResonance::calcLegendrePoly( const Double_t cosHel )
338 {
339  cosHel_ = cosHel;
340  return this->calcLegendrePoly();
341 }
342 
344 {
345  Double_t legPol = 1.0;
346 
347  if (resSpin_ == 1) {
348  // Calculate vector resonance Legendre polynomial
349  legPol = -2.0*cosHel_;
350  } else if (resSpin_ == 2) {
351  // Calculate tensor resonance Legendre polynomial
352  legPol = 4.0*(3.0*cosHel_*cosHel_ - 1.0)/3.0;
353  } else if (resSpin_ == 3) {
354  // Calculate spin 3 resonance Legendre polynomial
355  legPol = -8.0*(5.0*cosHel_*cosHel_*cosHel_ - 3.0*cosHel_)/5.0;
356  } else if (resSpin_ == 4) {
357  // Calculate spin 4 resonance Legendre polynomial
358  legPol = 16.0*(35.0*cosHel_*cosHel_*cosHel_*cosHel_ - 30.0*cosHel_*cosHel_ + 3.0)/35.0;
359  } else if (resSpin_ == 5) {
360  // Calculate spin 5 resonance Legendre polynomial
361  legPol = -32.0*(63.0*cosHel_*cosHel_*cosHel_*cosHel_*cosHel_ - 70.0*cosHel_*cosHel_*cosHel_ + 15.0*cosHel_)/63.0;
362  } else if (resSpin_ > 5) {
363  std::cerr << "WARNING in LauAbsResonance::calcLegendrePoly : Legendre polynomials not (yet) implemented for spin > 5" << std::endl;
364  }
365 
366  return legPol;
367 }
368 
369 void LauAbsResonance::changeResonance(const Double_t newMass, const Double_t newWidth, const Int_t newSpin)
370 {
371  if (newMass > 0.0) {
372  resMass_->valueAndRange(newMass,0.0,3.0*newMass);
373  resMass_->initValue(newMass);
374  resMass_->genValue(newMass);
375  std::cout << "INFO in LauAbsResonance::changeResonance : Setting mass to " << resMass_->value() << std::endl;
376  }
377  if (newWidth > 0.0) {
378  resWidth_->valueAndRange(newWidth,0.0,3.0*newWidth);
379  resWidth_->initValue(newWidth);
380  resWidth_->genValue(newWidth);
381  std::cout << "INFO in LauAbsResonance::changeResonance : Setting width to " << resWidth_->value() << std::endl;
382  }
383  if (newSpin > -1) {
384  resSpin_ = newSpin;
385  std::cout << "INFO in LauAbsResonance::changeResonance : Setting spin to " << resSpin_ << std::endl;
386  }
387 }
388 
389 void LauAbsResonance::changeBWBarrierRadii(const Double_t resRadius, const Double_t parRadius)
390 {
391  if ( resRadius >= 0.0 && resBWFactor_ != 0 ) {
392  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
393  resBWRadius->value(resRadius);
394  resBWRadius->initValue(resRadius);
395  resBWRadius->genValue(resRadius);
396  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting resonance factor radius to " << resBWRadius->value() << std::endl;
397  }
398  if ( parRadius >= 0.0 && parBWFactor_ != 0 ) {
399  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
400  parBWRadius->value(parRadius);
401  parBWRadius->initValue(parRadius);
402  parBWRadius->genValue(parRadius);
403  std::cout << "INFO in LauAbsResonance::changeBWBarrierRadii : Setting parent factor radius to " << parBWRadius->value() << std::endl;
404  }
405 }
406 
407 void LauAbsResonance::setResonanceParameter(const TString& name, const Double_t value)
408 {
409  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
410  std::cerr << "WARNING in LauAbsResonance::setResonanceParameter : Unable to set parameter \"" << name << "\" to value: " << value << "." << std::endl;
411 }
412 
414 {
415  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
416  std::cerr << "WARNING in LauAbsResonance::floatResonanceParameter : Unable to release parameter \"" << name << "\"." << std::endl;
417 }
418 
420 {
421  //This function should always be overwritten if needed in classes inheriting from LauAbsResonance.
422  std::cerr << "WARNING in LauAbsResonance::getResonanceParameter : Unable to get parameter \"" << name << "\"." << std::endl;
423  return 0;
424 }
425 
427 {
428  if ( param == 0 ) {
429  return;
430  }
431 
432  if ( param->clone() ) {
433  resParameters_.push_back( param->parent() );
434  } else {
435  resParameters_.push_back( param );
436  }
437 }
438 
439 void LauAbsResonance::fixBarrierRadii(const Bool_t fixResRad, const Bool_t fixParRad)
440 {
441  if ( resBWFactor_ == 0 ) {
442  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : resonance barrier factor not present, cannot fix/float it" << std::endl;
443  return;
444  }
445 
446  if ( parBWFactor_ == 0 ) {
447  std::cerr << "WARNING in LauAbsResonance::fixBarrierRadii : parent barrier factor not present, cannot fix/float it" << std::endl;
448  return;
449  }
450 
451  LauParameter* resBWRadius = resBWFactor_->getRadiusParameter();
452  resBWRadius->fixed(fixResRad);
453 
454  LauParameter* parBWRadius = parBWFactor_->getRadiusParameter();
455  parBWRadius->fixed(fixParRad);
456 }
457 
459 {
460  if ( resBWFactor_ == 0 ) {
461  std::cerr << "WARNING in LauAbsResonance::fixResRadius : resonance barrier factor not present" << std::endl;
462  return kTRUE;
463  }
464 
466  return bwRadius->fixed();
467 }
468 
470 {
471  if ( parBWFactor_ == 0 ) {
472  std::cerr << "WARNING in LauAbsResonance::fixParRadius : parent barrier factor not present" << std::endl;
473  return kTRUE;
474  }
475 
477  return bwRadius->fixed();
478 }
479 
481 {
482  if ( resBWFactor_ == 0 ) {
483  std::cerr << "WARNING in LauAbsResonance::getResRadius : resonance barrier factor not present" << std::endl;
484  return -1.0;
485  }
486 
488  return bwRadius->unblindValue();
489 }
490 
492 {
493  if ( parBWFactor_ == 0 ) {
494  std::cerr << "WARNING in LauAbsResonance::getParRadius : parent barrier factor not present" << std::endl;
495  return -1.0;
496  }
497 
499  return bwRadius->unblindValue();
500 }
501 
503 {
504  // Get the parent mass
505  Double_t mass(LauConstants::mB);
506 
507  if (daughters_) {
508  mass = daughters_->getMassParent();
509  }
510 
511  return mass;
512 }
513 
515 {
516  // Get the daughter mass
517  Double_t mass(LauConstants::mPi);
518 
519  if (daughters_) {
520  if (resPairAmpInt_ == 1) {
521  mass = daughters_->getMassDaug2();
522  } else if (resPairAmpInt_ == 2) {
523  mass = daughters_->getMassDaug1();
524  } else if (resPairAmpInt_ == 3) {
525  mass = daughters_->getMassDaug1();
526  }
527  }
528 
529  return mass;
530 }
531 
533 {
534  // Get the daughter mass
535  Double_t mass(LauConstants::mPi);
536 
537  if (daughters_) {
538  if (resPairAmpInt_ == 1) {
539  mass = daughters_->getMassDaug3();
540  } else if (resPairAmpInt_ == 2) {
541  mass = daughters_->getMassDaug3();
542  } else if (resPairAmpInt_ == 3) {
543  mass = daughters_->getMassDaug2();
544  }
545  }
546 
547  return mass;
548 }
549 
551 {
552  // Get the bachelor mass
553  Double_t mass(LauConstants::mPi);
554 
555  if (daughters_) {
556  if (resPairAmpInt_ == 1) {
557  mass = daughters_->getMassDaug1();
558  } else if (resPairAmpInt_ == 2) {
559  mass = daughters_->getMassDaug2();
560  } else if (resPairAmpInt_ == 3) {
561  mass = daughters_->getMassDaug3();
562  }
563  }
564 
565  return mass;
566 }
567 
569 {
570  // Get the parent charge
571  Int_t charge(0);
572 
573  if (daughters_) {
574  charge = daughters_->getChargeParent();
575  }
576 
577  return charge;
578 }
579 
581 {
582  // Get the daughter charge
583  Int_t charge(0);
584 
585  if (daughters_) {
586  if (resPairAmpInt_ == 1) {
587  charge = daughters_->getChargeDaug2();
588  } else if (resPairAmpInt_ == 2) {
589  charge = daughters_->getChargeDaug1();
590  } else if (resPairAmpInt_ == 3) {
591  charge = daughters_->getChargeDaug1();
592  }
593  }
594 
595  return charge;
596 }
597 
599 {
600  // Get the daughter charge
601  Int_t charge(0);
602 
603  if (daughters_) {
604  if (resPairAmpInt_ == 1) {
605  charge = daughters_->getChargeDaug3();
606  } else if (resPairAmpInt_ == 2) {
607  charge = daughters_->getChargeDaug3();
608  } else if (resPairAmpInt_ == 3) {
609  charge = daughters_->getChargeDaug2();
610  }
611  }
612 
613  return charge;
614 }
615 
617 {
618  // Get the bachelor charge
619  Int_t charge(0);
620 
621  if (daughters_) {
622  if (resPairAmpInt_ == 1) {
623  charge = daughters_->getChargeDaug1();
624  } else if (resPairAmpInt_ == 2) {
625  charge = daughters_->getChargeDaug2();
626  } else if (resPairAmpInt_ == 3) {
627  charge = daughters_->getChargeDaug3();
628  }
629  }
630 
631  return charge;
632 }
633 
635 {
636  // Get the parent name
637  TString name("");
638 
639  if (daughters_) {
640  name = daughters_->getNameParent();
641  }
642 
643  return name;
644 }
645 
647 {
648  // Get the daughter name
649  TString name("");
650 
651  if (daughters_) {
652  if (resPairAmpInt_ == 1) {
653  name = daughters_->getNameDaug2();
654  } else if (resPairAmpInt_ == 2) {
655  name = daughters_->getNameDaug1();
656  } else if (resPairAmpInt_ == 3) {
657  name = daughters_->getNameDaug1();
658  }
659  }
660 
661  return name;
662 }
663 
665 {
666  // Get the daughter name
667  TString name("");
668 
669  if (daughters_) {
670  if (resPairAmpInt_ == 1) {
671  name = daughters_->getNameDaug3();
672  } else if (resPairAmpInt_ == 2) {
673  name = daughters_->getNameDaug3();
674  } else if (resPairAmpInt_ == 3) {
675  name = daughters_->getNameDaug2();
676  }
677  }
678 
679  return name;
680 }
681 
683 {
684  // Get the bachelor name
685  TString name("");
686 
687  if (daughters_) {
688  if (resPairAmpInt_ == 1) {
689  name = daughters_->getNameDaug1();
690  } else if (resPairAmpInt_ == 2) {
691  name = daughters_->getNameDaug2();
692  } else if (resPairAmpInt_ == 3) {
693  name = daughters_->getNameDaug3();
694  }
695  }
696 
697  return name;
698 }
699 
Double_t getMassParent() const
Get mass of the parent particle.
TString getNameDaug1() const
Get name of the first daughter particle.
Double_t mass_
Invariant mass.
Int_t getChargeDaug3() const
Get charge of the third daughter particle.
Double_t calcCovSpinFactor(const Double_t pProd)
Calculate the amplitude spin term using the covariant tensor formalism.
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 calcCovFactor(const Double_t erm)
Calculate the spin-dependent covariant factor.
Double_t calcLegendrePoly() const
Calculate the Legendre polynomial for the spin factor.
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 getcov12() const
Get covariant factor in 12 axis.
Double_t getm23() const
Get the m23 invariant mass.
Double_t massDaug1_
Daughter 1 mass.
Class for defining the properties of a resonant particle.
Double_t getcov23() const
Get covariant factor in 23 axis.
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.
Double_t covFactor_
Covariant factor (full spin-dependent expression)
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 calcZemachSpinFactor(const Double_t pProd) const
Calculate the amplitude spin term using the Zemach tensor formalism.
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.
Double_t getcov13() const
Get covariant factor in 13 axis.
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.
Double_t erm_
Covariant factor.
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 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
LauSpinType getSpinType() const
Get the spin type.
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.
Double_t cosHel_
Helicity angle cosine.
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.