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