laura is hosted by Hepforge, IPPP Durham
Laura++  v3r1
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauKinematics.cc
Go to the documentation of this file.
1 
2 // Copyright University of Warwick 2004 - 2013.
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 "TF2.h"
18 #include "TMath.h"
19 #include "TRandom.h"
20 
21 #include "LauConstants.hh"
22 #include "LauKinematics.hh"
23 #include "LauRandom.hh"
24 
26 
27 
28 LauKinematics::LauKinematics(Double_t m1, Double_t m2, Double_t m3, Double_t mParent, Bool_t calcSquareDPCoords) :
29  m1_(m1), m2_(m2), m3_(m3), mParent_(mParent),
30  m1Sq_(m1*m1), m2Sq_(m2*m2), m3Sq_(m3*m3), mParentSq_(mParent*mParent),
31  mDTot_(m1 + m2 + m3),
32  massInt_(mParent - (m1+m2+m3)),
33  mSqDTot_(m1*m1 + m2*m2 + m3*m3),
34  m12_(0.0), m23_(0.0), m13_(0.0),
35  m12Sq_(0.0), m23Sq_(0.0), m13Sq_(0.0),
36  c12_(0.0), c23_(0.0), c13_(0.0),
37  mPrime_(0.0), thetaPrime_(0.0),
38  qi_(0.0), qk_(0.0),
39  p1_12_(0.0), p3_12_(0.0),
40  p2_23_(0.0), p1_23_(0.0),
41  p1_13_(0.0), p2_13_(0.0),
42  p1_Parent_(0.0), p2_Parent_(0.0), p3_Parent_(0.0),
43  squareDP_(calcSquareDPCoords), warnings_(kTRUE)
44 {
45  // Constructor
46  mass_.clear(); mMin_.clear(); mMax_.clear(); mSqMin_.clear(); mSqMax_.clear();
47  mSq_.clear();
48  mSqDiff_.clear();
49  mDiff_.clear();
50 
51  mass_.push_back(m1_);
52  mass_.push_back(m2_);
53  mass_.push_back(m3_);
54 
55  mSq_.push_back(m1Sq_);
56  mSq_.push_back(m2Sq_);
57  mSq_.push_back(m3Sq_);
58 
59  // DP max and min kinematic boundary for circumscribed box
60  // (see figure in PDG book).
61  for (Int_t i = 0; i < 3; i++) {
62  mMin_.push_back(mDTot_ - mass_[i]);
63  mMax_.push_back(mParent_ - mass_[i]);
64  mSqMin_.push_back(mMin_[i]*mMin_[i]);
65  mSqMax_.push_back(mMax_[i]*mMax_[i]);
66  mSqDiff_.push_back(mSqMax_[i] - mSqMin_[i]);
67  mDiff_.push_back(mMax_[i] - mMin_[i]);
68  }
69 
70  if (this->squareDP()) {
71  std::cout<<"INFO in LauKinematics::LauKinematics : squareDP = kTRUE"<<std::endl;
72  } else {
73  std::cout<<"INFO in LauKinematics::LauKinematics : squareDP = kFALSE"<<std::endl;
74  }
75 }
76 
78 {
79  // Destructor
80 }
81 
82 void LauKinematics::updateKinematics(Double_t m13Sq, Double_t m23Sq)
83 {
84  // Sets the internal private data members
85  // m13Sq and m23Sq, as well as m13 and m23.
86  // Also calculate m12Sq and m12, given mParent defined in the constructor.
87  // Lastly, calculate the helicity cosines.
88 
89  // Update the 3 mass-squares
90  this->updateMassSquares(m13Sq, m23Sq);
91 
92  // Now update the helicity cosines
93  this->calcHelicities();
94 
95  // Also calculate the m', theta' variables
96  if (squareDP_) {this->calcSqDPVars();}
97 }
98 
99 void LauKinematics::updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
100 {
101  // From square DP co-ordinates, calculate remaining kinematic variables
102  this->updateSqDPMassSquares(mPrime, thetaPrime);
103 
104  // Finally calculate the helicities and track momenta
105  this->calcHelicities();
106 }
107 
109 {
110  // For given m_12 and cos(theta_12) values, calculate m' and theta' for the square Dalitz plot
111  Double_t value = (2.0*(m12_ - mMin_[2])/mDiff_[2]) - 1.0;
112  mPrime_ = LauConstants::invPi*TMath::ACos(value);
113  thetaPrime_ = LauConstants::invPi*TMath::ACos(c12_);
114  // Sometimes events are assigned exactly thetaPrime = 0.0 or 1.0
115  // which gives problems with efficiency and other histograms
116  if (thetaPrime_ == 0.0)
117  {
118  thetaPrime_ += 1.0e-10;
119  } else if (thetaPrime_ == 1.0)
120  {
121  thetaPrime_ -= 1.0e-10;
122  }
123 }
124 
126 {
127  return this->calcSqDPJacobian(mPrime_,thetaPrime_);
128 }
129 
130 Double_t LauKinematics::calcSqDPJacobian(const Double_t mPrime, const Double_t thPrime) const
131 {
132  // Calculate the Jacobian for the transformation
133  // m23^2, m13^2 -> m', theta' (square DP)
134  const Double_t m12 = 0.5*mDiff_[2]*(1.0 + TMath::Cos(LauConstants::pi*mPrime)) + mMin_[2];
135  const Double_t m12Sq = m12*m12;
136 
137  const Double_t e1Cms12 = (m12Sq - m2Sq_ + m1Sq_)/(2.0*m12);
138  const Double_t e3Cms12 = (mParentSq_ - m12Sq - m3Sq_)/(2.0*m12);
139 
140  const Double_t p1Cms12 = this->pCalc(e1Cms12, m1Sq_);
141  const Double_t p3Cms12 = this->pCalc(e3Cms12, m3Sq_);
142 
143  const Double_t deriv1 = LauConstants::piBy2*mDiff_[2]*TMath::Sin(LauConstants::pi*mPrime);
144  const Double_t deriv2 = LauConstants::pi*TMath::Sin(LauConstants::pi*thPrime);
145 
146  const Double_t jacobian = 4.0*p1Cms12*p3Cms12*m12*deriv1*deriv2;
147 
148  return jacobian;
149 }
150 
151 void LauKinematics::updateMassSquares(Double_t m13Sq, Double_t m23Sq)
152 {
153  m13Sq_ = m13Sq;
154  if (m13Sq_ > 0.0) {
155  m13_ = TMath::Sqrt(m13Sq_);
156  } else {
157  m13_ = 0.0;
158  }
159 
160  m23Sq_ = m23Sq;
161  if (m23Sq_ > 0.0) {
162  m23_ = TMath::Sqrt(m23Sq_);
163  } else {
164  m23_ = 0.0;
165  }
166 
167  // Now calculate m12Sq and m12.
168  this->calcm12Sq();
169 
170  // Calculate momenta of tracks in parent (B, D etc.) rest-frame
171  this->calcParentFrameMomenta();
172 }
173 
174 void LauKinematics::updateSqDPMassSquares(Double_t mPrime, Double_t thetaPrime)
175 {
176  // From square DP co-ordinates, calculate only the mass-squares
177 
178  // First set the square-DP variables
179  mPrime_ = mPrime; thetaPrime_ = thetaPrime;
180 
181  // Next calculate m12 and c12
182  Double_t m12 = 0.5*mDiff_[2]*(1.0 + TMath::Cos(LauConstants::pi*mPrime)) + mMin_[2];
183  Double_t c12 = TMath::Cos(LauConstants::pi*thetaPrime);
184 
185  // From these calculate m13Sq and m23Sq
186  this->updateMassSq_m12(m12, c12);
187 
188  // Calculate momenta of tracks in parent (B, D etc.) rest-frame
189  this->calcParentFrameMomenta();
190 }
191 
193 {
194  // Calculate m12Sq from m13Sq and m23Sq.
196 
197  // If m12Sq is too low, return lower limit,
198  // and modify m13Sq accordingly.
199  if (m12Sq_ < mSqMin_[2]) {
200  m12Sq_ = mSqMin_[2] + 1.0e-3;
202  }
203 
204  if (m12Sq_ > 0.0) {
205  m12_ = TMath::Sqrt(m12Sq_);
206  } else {
207  m12_ = 0.0;
208  }
209 }
210 
212 {
213  Double_t e1 = (mParentSq_ + m1Sq_ - m23Sq_) / (2.0*mParent_);
214  Double_t e2 = (mParentSq_ + m2Sq_ - m13Sq_) / (2.0*mParent_);
215  Double_t e3 = (mParentSq_ + m3Sq_ - m12Sq_) / (2.0*mParent_);
216 
217  p1_Parent_ = this->pCalc(e1, m1Sq_);
218  p2_Parent_ = this->pCalc(e2, m2Sq_);
219  p3_Parent_ = this->pCalc(e3, m3Sq_);
220 }
221 
223 {
224  // Calculate helicity angle cosines, given m12Sq, m13Sq and m23Sq.
225  // The notation is confusing for the helicity angle cosines:
226  // cij is helicity angle for the pair which is made from tracks i and j.
227  // It is the angle between tracks i and k in the ij rest frame
228  // Make sure that setMassSqVars is run first.
229  Int_t zero(0), one(1), two(2);
230 
231  c12_ = cFromM(m12Sq_, m13Sq_, m12_, zero, one, two);
232  p1_12_ = qi_; p3_12_ = qk_; // i, j = 12 (rest frame), k = 3
233  c23_ = cFromM(m23Sq_, m12Sq_, m23_, one, two, zero);
234  p2_23_ = qi_; p1_23_ = qk_; // i, j = 23 (rest frame), k = 1
235  c13_ = cFromM(m13Sq_, m23Sq_, m13_, two, zero, one);
236  p1_13_ = qi_; p2_13_ = qk_; // i, j = 31 (rest frame), k = 2
237 
238 }
239 
240 Double_t LauKinematics::cFromM(Double_t mijSq, Double_t mikSq, Double_t mij, Int_t i, Int_t j, Int_t k) const
241 {
242  // Routine to calculate the cos(helicity) variables from the masses of the particles.
243  // cij is helicity angle for the pair which is made from tracks i and j.
244  // It is the angle between tracks i and k in the ij rest frame.
245  Double_t EiCmsij = (mijSq - mSq_[j] + mSq_[i])/(2.0*mij);
246  Double_t EkCmsij = (mParentSq_ - mijSq - mSq_[k])/(2.0*mij);
247 
248  if (EiCmsij < mass_[i]) {
249  if (warnings_) {
250  std::cerr<<"WARNING in LauKinematics::cFromM : EiCmsij = "<<EiCmsij<<" too small < mass_["<<i<<"] = "<<mass_[i]<<" in cFromM, i, j, k = "<<i<<", "<<j<<", "<<k<<std::endl;
251  std::cerr<<" : mijSq = "<<mijSq<<"; mij = "<<mij<<"; mSq_["<<j<<"] = "<<mSq_[j]<<"; mSq_["<<i<<"] = "<<mSq_[i]<<std::endl;
252  }
253  return 0.0;
254  }
255  if (EkCmsij < mass_[k]) {
256  if (warnings_) {
257  std::cerr<<"WARNING in LauKinematics::cFromM : EkCmsij = "<<EkCmsij<<" too small < mass_["<<k<<"] = "<<mass_[k]<<" in cFromM, i, j, k = "<<i<<", "<<j<<", "<<k<<std::endl;
258  std::cerr<<" : mijSq = "<<mijSq<<"; mij = "<<mij<<"; mSq_["<<j<<"] = "<<mSq_[j]<<"; mSq_["<<i<<"] = "<<mSq_[i]<<std::endl;
259  }
260  return 0.0;
261  }
262 
263  // Find track i and k momenta in ij rest frame
264  qi_ = this->pCalc(EiCmsij, mSq_[i]);
265  qk_ = this->pCalc(EkCmsij, mSq_[k]);
266 
267  // Find ij helicity angle
268  Double_t cosHel = -(mikSq - mSq_[i] - mSq_[k] - 2.0*EiCmsij*EkCmsij)/(2.0*qi_*qk_);
269 
270  if (cosHel > 1.0) {
271  cosHel = 1.0;
272  } else if (cosHel < -1.0) {
273  cosHel = -1.0;
274  }
275 
276  if (i == 1) {cosHel *= -1.0;}
277 
278  return cosHel;
279 
280 }
281 
282 Double_t LauKinematics::mFromC(Double_t mijSq, Double_t cij, Double_t mij, Int_t i, Int_t j, Int_t k) const
283 {
284  // Returns the mass-squared for a pair of particles, i,j, given their
285  // invariant mass (squared) and the helicity angle.
286  // cij is helicity angle for the pair which is made from tracks i and j.
287  // It is the angle between tracks i and k in the ij rest frame.
288 
289  Double_t EiCmsij = (mijSq - mSq_[j] + mSq_[i])/(2.0*mij);
290  Double_t EkCmsij = (mParentSq_ - mijSq - mSq_[k])/(2.0*mij);
291 
292  if (TMath::Abs(EiCmsij - mass_[i]) > 1e-6 && EiCmsij < mass_[i]) {
293  if (warnings_) {
294  std::cerr<<"WARNING in LauKinematics::mFromC : EiCmsij = "<<EiCmsij<<" < "<<mass_[i]<<" in mFromC, i, j, k = "<<i<<", "<<j<<", "<<k<<std::endl;
295  }
296  return 0.0;
297  }
298  if (TMath::Abs(EkCmsij - mass_[k]) > 1e-6 && EkCmsij < mass_[k]) {
299  if (warnings_) {
300  std::cerr<<"WARNING in LauKinematics::mFromC : EkCmsij = "<<EkCmsij<<" < "<<mass_[k]<<" in mFromC, i, j, k = "<<i<<", "<<j<<", "<<k<<std::endl;
301  }
302  return 0.0;
303  }
304 
305  // Find track i and k momenta in ij rest fram
306  Double_t qi = this->pCalc(EiCmsij, mSq_[i]);
307  Double_t qk = this->pCalc(EkCmsij, mSq_[k]);
308 
309  // Find mikSq
310  Double_t massSq = mSq_[i] + mSq_[k] + 2.0*EiCmsij*EkCmsij - 2.0*qi*qk*cij;
311 
312  if (massSq < mSqMin_[j]) {
313  if (warnings_) {
314  std::cerr<<"WARNING in LauKinematics::mFromC : mFromC below bound: i, j, k = "<<i<<", "<<j<<", "<<k<<std::endl;
315  }
316  massSq = mSqMin_[j];
317  }
318 
319  return massSq;
320 
321 }
322 
323 void LauKinematics::genFlatPhaseSpace(Double_t& m13Sq, Double_t& m23Sq) const
324 {
325  // Routine to generate flat phase-space events.
326  // DP max kinematic boundaries in circumscribed box
327  // See DP figure in PDG book.
328  // m13max=mbrec-mass(2)
329  // m13sqmax=m13max*m13max
330  // m23max=mbrec-mass(1)
331  // m23sqmax=m23max*m23max
332 
333  // Generate m13Sq and m23Sq flat within DP overall rectangular box
334  // Loop if the generated point is not within the DP kinematic boundary
335 
336  do {
337  m13Sq = mSqMin_[1] + LauRandom::randomFun()->Rndm()*mSqDiff_[1];
338  m23Sq = mSqMin_[0] + LauRandom::randomFun()->Rndm()*mSqDiff_[0];
339 
340  } while ( ! this->withinDPLimits( m13Sq, m23Sq ) );
341 }
342 
343 void LauKinematics::genFlatSqDP(Double_t& mPrime, Double_t& thetaPrime) const
344 {
345  // Generate random event in the square Dalitz plot
346  mPrime = LauRandom::randomFun()->Rndm();
347  thetaPrime = LauRandom::randomFun()->Rndm();
348 }
349 
350 Bool_t LauKinematics::withinDPLimits(Double_t m13Sq, Double_t m23Sq) const
351 {
352  // Find out whether the point (m13Sq,m23Sq) is within the limits of the
353  // Dalitz plot. The limits are specified by the invariant masses
354  // of the parent (e.g. B) and its three daughters that were
355  // defined in the constructor of this class. Here
356  // m_13Sq = square of invariant mass of daughters 1 and 3
357  // m_23Sq = square of invariant mass of daughters 2 and 3.
358 
359  Bool_t withinDP = kFALSE;
360 
361  // First check that m13Sq is within its absolute min and max
362  if (!((m13Sq > mSqMin_[1]) && (m13Sq < mSqMax_[1]))) {
363  return kFALSE;
364  }
365 
366  // Now for the given value of m13Sq calculate the local min and max of m23Sq
367  Double_t m13 = TMath::Sqrt(m13Sq);
368 
369  Double_t e3Cms13 = (m13Sq - m1Sq_ + m3Sq_)/(2.0*m13);
370  Double_t p3Cms13 = this->pCalc(e3Cms13, m3Sq_);
371 
372  Double_t e2Cms13 = (mParentSq_ - m13Sq - m2Sq_)/(2.0*m13);
373  Double_t p2Cms13 = this->pCalc(e2Cms13, m2Sq_);
374 
375  Double_t term = 2.0*e2Cms13*e3Cms13 + m2Sq_ + m3Sq_;
376 
377  Double_t m23SqLocMin = term - 2.0*p2Cms13*p3Cms13;
378  Double_t m23SqLocMax = term + 2.0*p2Cms13*p3Cms13;
379 
380  // Check whether the given value of m23Sq satisfies these bounds
381  if (m23Sq > m23SqLocMin && m23Sq < m23SqLocMax) {
382  withinDP = kTRUE;
383  }
384 
385  return withinDP;
386 }
387 
388 Bool_t LauKinematics::withinDPLimits2(Double_t m13Sq, Double_t m23Sq) const
389 {
390  // Same as withinDPLimits, but this time testing whether the m13Sq
391  // variable is within the kinematic boundary for the given m23Sq value
392 
393  Bool_t withinDP = kFALSE;
394 
395  // First check that m13Sq is within its absolute min and max
396  if (!((m23Sq > mSqMin_[0]) && (m23Sq < mSqMax_[0]))) {
397  return kFALSE;
398  }
399 
400  // Now for the given value of m13Sq calculate the local min and max of m23Sq
401  Double_t m23 = TMath::Sqrt(m23Sq);
402 
403  Double_t e3Cms23 = (m23Sq - m2Sq_ + m3Sq_)/(2.0*m23);
404  Double_t p3Cms23 = this->pCalc(e3Cms23, m3Sq_);
405 
406  Double_t e1Cms23 = (mParentSq_ - m23Sq - m1Sq_)/(2.0*m23);
407  Double_t p1Cms23 = this->pCalc(e1Cms23, m1Sq_);
408 
409  Double_t term = 2.0*e3Cms23*e1Cms23 + m1Sq_ + m3Sq_;
410 
411  Double_t m13SqLocMin = term - 2.0*p3Cms23*p1Cms23;
412  Double_t m13SqLocMax = term + 2.0*p3Cms23*p1Cms23;
413 
414  // Check whether the given value of m23Sq satisfies these bounds
415  if (m13Sq > m13SqLocMin && m13Sq < m13SqLocMax) {
416  withinDP = kTRUE;
417  }
418 
419  return withinDP;
420 }
421 
422 Bool_t LauKinematics::withinSqDPLimits(Double_t mPrime, Double_t thetaPrime) const
423 {
424  // Check whether m' and theta' are between 0 and 1
425  Bool_t withinDP(kFALSE);
426  if (mPrime > 0.0 && mPrime < 1.0 &&
427  thetaPrime > 0.0 && thetaPrime < 1.0) {
428  withinDP = kTRUE;
429  }
430 
431  return withinDP;
432 }
433 
434 Double_t LauKinematics::calcThirdMassSq(Double_t firstMassSq, Double_t secondMassSq) const
435 {
436  // Calculate one massSq from the other two
437  return mParentSq_ + mSqDTot_ - firstMassSq - secondMassSq;
438 }
439 
441 {
442  return this->distanceFromDPCentre(m13Sq_,m23Sq_);
443 }
444 
445 Double_t LauKinematics::distanceFromDPCentre(Double_t m13Sq, Double_t m23Sq) const
446 {
447  // DP centre is defined as the point where m12=m13=m23
448 
449  // First find the m^2_ij value for the centre
450  Double_t centreMijSq = (mParentSq_ + m1Sq_ + m2Sq_ + m3Sq_)/3.0;
451 
452  // Then find the difference between this and the two provided co-ords
453  Double_t diff13 = m13Sq - centreMijSq;
454  Double_t diff23 = m23Sq - centreMijSq;
455 
456  // Calculate the total distance
457  Double_t distance = TMath::Sqrt(diff13*diff13 + diff23*diff23);
458  return distance;
459 }
460 
461 Double_t LauKinematics::pCalc(Double_t energy, Double_t massSq) const
462 {
463  // Routine to calculate the momentum of a particle, given its energy and
464  // invariant mass (squared).
465  Double_t arg = energy*energy - massSq;
466 
467  if (arg < 0.0) {
468  //if (warnings_) {
469  //std::cerr<<"WARNING in LauKinematics::pCalc : arg < 0.0: "<<arg<<" for e = "<<energy<<", mSq = "<<massSq<<std::endl;
470  //}
471  arg = 0.0;
472  }
473  Double_t pCalcVal = TMath::Sqrt(arg);
474  return pCalcVal;
475 
476 }
477 
479 {
480  // Flips the DP variables m13^2 <-> m23^2.
481  // Used in the case of symmetrical Dalitz plots (i.e. when two of the
482  // daughter tracks are the same type) within the
483  // LauIsobarDynamics::resAmp function.
485 }
486 
488 {
489  // Cyclically rotates the DP variables (m12 -> m23, m23 -> m13, m13 -> m12)
490  // Used in the case of fully symmetric Dalitz plots (i.e. when all
491  // three of the daughter tracks are the same type) within the
492  // LauIsobarDynamics::resAmp function.
494 }
495 
496 void LauKinematics::updateMassSq_m23(Double_t m23, Double_t c23)
497 {
498  // Update the variables m12Sq_ and m13Sq_ given m23 and c23.
499  m23_ = m23; m23Sq_ = m23*m23; c23_ = c23;
500 
501  Int_t zero(0), one(1), two(2);
502  m13Sq_ = this->mFromC(m23Sq_, -c23_, m23_, two, one, zero);
504 
505  m13_ = TMath::Sqrt(m13Sq_);
506  m12_ = TMath::Sqrt(m12Sq_);
507  //c12_ = cFromM(m12Sq_, m13Sq_, m12_, zero, one, two);
508  //c13_ = cFromM(m13Sq_, m23Sq_, m13_, two, zero, one);
509 }
510 
511 void LauKinematics::updateMassSq_m13(Double_t m13, Double_t c13)
512 {
513  // Update the variables m12Sq_ and m23Sq_ given m13 and c13.
514  m13_ = m13; m13Sq_ = m13*m13; c13_ = c13;
515 
516  Int_t zero(0), one(1), two(2);
517  m23Sq_ = this->mFromC(m13Sq_, c13_, m13_, two, zero, one);
519 
520  m23_ = TMath::Sqrt(m23Sq_);
521  m12_ = TMath::Sqrt(m12Sq_);
522  //c23_ = cFromM(m23Sq_, m12Sq_, m23_, one, two, zero);
523  //c12_ = cFromM(m12Sq_, m13Sq_, m12_, zero, one, two);
524 }
525 
526 void LauKinematics::updateMassSq_m12(Double_t m12, Double_t c12)
527 {
528  // Update the variables m23Sq_ and m13Sq_ given m12 and c12.
529  m12_ = m12; m12Sq_ = m12*m12; c12_ = c12;
530 
531  Int_t zero(0), one(1), two(2);
532  m13Sq_ = this->mFromC(m12Sq_, c12_, m12_, zero, one, two);
534  m13_ = TMath::Sqrt(m13Sq_);
535  m23_ = TMath::Sqrt(m23Sq_);
536  //c23_ = cFromM(m23Sq_, m12Sq_, m23_, one, two, zero);
537  //c13_ = cFromM(m13Sq_, m23Sq_, m13_, two, zero, one);
538 }
539 
540 Double_t LauKinematics::genm13Sq() const
541 {
542  Double_t m13Sq = mSqMin_[1] + LauRandom::randomFun()->Rndm()*mSqDiff_[1];
543  return m13Sq;
544 }
545 
546 Double_t LauKinematics::genm23Sq() const
547 {
548  Double_t m23Sq = mSqMin_[0] + LauRandom::randomFun()->Rndm()*mSqDiff_[0];
549  return m23Sq;
550 }
551 
552 Double_t LauKinematics::genm12Sq() const
553 {
554  Double_t m12Sq = mSqMin_[2] + LauRandom::randomFun()->Rndm()*mSqDiff_[2];
555  return m12Sq;
556 }
557 
558 
559 void LauKinematics::drawDPContour(Int_t orientation, Int_t nbins)
560 {
561  // orientation -
562  // 1323 : x = m13, y = m23
563  // etc.
564 
565  Double_t m1 = this->getm1();
566  Double_t m2 = this->getm2();
567  Double_t m3 = this->getm3();
568  Double_t mParent = this->getmParent();
569 
570  Double_t m13SqMin = this->getm13SqMin();
571  Double_t m23SqMin = this->getm23SqMin();
572  Double_t m12SqMin = this->getm12SqMin();
573  Double_t m13SqMax = this->getm13SqMax();
574  Double_t m23SqMax = this->getm23SqMax();
575  Double_t m12SqMax = this->getm12SqMax();
576 
577  Double_t xMin(0.0);
578  Double_t xMax(mParent*mParent);
579  Double_t yMin(0.0);
580  Double_t yMax(mParent*mParent);
581  if (orientation == 1323) {
582  xMin = m13SqMin-1.0; xMax = m13SqMax+1.0;
583  yMin = m23SqMin-1.0; yMax = m23SqMax+1.0;
584  } else if (orientation == 2313) {
585  xMin = m23SqMin-1.0; xMax = m23SqMax+1.0;
586  yMin = m13SqMin-1.0; yMax = m13SqMax+1.0;
587  } else if (orientation == 1213) {
588  xMin = m12SqMin-1.0; xMax = m12SqMax+1.0;
589  yMin = m13SqMin-1.0; yMax = m13SqMax+1.0;
590  } else if (orientation == 1312) {
591  xMin = m13SqMin-1.0; xMax = m13SqMax+1.0;
592  yMin = m12SqMin-1.0; yMax = m12SqMax+1.0;
593  } else if (orientation == 1223) {
594  xMin = m12SqMin-1.0; xMax = m12SqMax+1.0;
595  yMin = m23SqMin-1.0; yMax = m23SqMax+1.0;
596  } else if (orientation == 2312) {
597  xMin = m23SqMin-1.0; xMax = m23SqMax+1.0;
598  yMin = m12SqMin-1.0; yMax = m12SqMax+1.0;
599  } else {
600  std::cerr<<"ERROR in LauKinematics::drawDPContour : Unrecognised orientation, not drawing contour."<<std::endl;
601  return;
602  }
603 
604  Int_t npar = 4;
605  TF2 * f2 = new TF2("contour", &dal, xMin, xMax, yMin, yMax, npar);
606 
607  // Set the parameters
608  f2->SetParameter(0,mParent);
609  if (orientation == 1323) {
610  f2->SetParameter(1,m1);
611  f2->SetParameter(2,m2);
612  f2->SetParameter(3,m3);
613  } else if (orientation == 2313) {
614  f2->SetParameter(1,m2);
615  f2->SetParameter(2,m1);
616  f2->SetParameter(3,m3);
617  } else if (orientation == 1213) {
618  f2->SetParameter(1,m2);
619  f2->SetParameter(2,m3);
620  f2->SetParameter(3,m1);
621  } else if (orientation == 1312) {
622  f2->SetParameter(1,m3);
623  f2->SetParameter(2,m2);
624  f2->SetParameter(3,m1);
625  } else if (orientation == 1223) {
626  f2->SetParameter(1,m1);
627  f2->SetParameter(2,m3);
628  f2->SetParameter(3,m2);
629  } else if (orientation == 2312) {
630  f2->SetParameter(1,m3);
631  f2->SetParameter(2,m1);
632  f2->SetParameter(3,m2);
633  }
634 
635  // Set up the contour to be drawn when the value of the function == 1.0
636  Double_t b[]={1.0};
637  f2->SetContour(1,b);
638 
639  // Set the number of bins for the contour to be sampled over
640  f2->SetNpx(nbins);
641  f2->SetNpy(nbins);
642  // and the line style
643  f2->SetLineWidth(3);
644  f2->SetLineStyle(kSolid);
645 
646  // Draw the contour on top of the histo that should already have been drawn
647  f2->DrawCopy("same");
648 
649  delete f2;
650 }
651 
652 Double_t dal(Double_t* x, Double_t* par)
653 {
654  Double_t mParent = par[0];
655  Double_t mi = par[1];
656  Double_t mj = par[2];
657  Double_t mk = par[3];
658 
659  Double_t mikSq=x[0];
660  Double_t mjkSq=x[1];
661  Double_t mik = TMath::Sqrt(mikSq);
662  Double_t mjk = TMath::Sqrt(mjkSq);
663 
664  Double_t ejcmsik = (mParent*mParent-mj*mj-mik*mik)/(2.0*mik);
665  Double_t ekcmsik = (mik*mik+mk*mk-mi*mi)/(2.0*mik);
666  if (ekcmsik<mk || ejcmsik<mj) return 2.0;
667 
668  Double_t pj = TMath::Sqrt(ejcmsik*ejcmsik-mj*mj);
669  Double_t pk = TMath::Sqrt(ekcmsik*ekcmsik-mk*mk);
670  Double_t coshelik = (mjk*mjk - mk*mk - mj*mj - 2.0*ejcmsik*ekcmsik)/(2.0*pj*pk);
671 
672  Double_t coshelikSq = coshelik*coshelik;
673  return coshelikSq;
674 }
675 
Double_t p1_23_
Momentum of track 1 in 2-3 rest frame.
Double_t getm13SqMin() const
Get the m13Sq minimum, (m1 + m3)^2.
const Double_t mParentSq_
Mass of parent particle squared.
Double_t calcThirdMassSq(Double_t firstMassSq, Double_t secondMassSq) const
Calculate the third invariant mass square from the two provided (e.g. mjkSq from mijSq and mikSq) ...
TRandom * randomFun()
Access the singleton random number generator with a particular seed.
Definition: LauRandom.cc:20
void updateMassSq_m13(Double_t m13, Double_t c13)
Update the variables m12Sq_ and m23Sq_ given the invariant mass m13 and the cosine of the helicity an...
Double_t getm23SqMin() const
Get the m23Sq minimum, (m2 + m3)^2.
Double_t qi_
Momentum q of particle i.
void updateSqDPKinematics(Double_t mPrime, Double_t thetaPrime)
Update all kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
std::vector< Double_t > mSqMin_
Vector of the minimum mijSq values.
ClassImp(LauAbsCoeffSet)
Bool_t withinDPLimits2(Double_t m13Sq, Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot (alternat...
Double_t m23_
Invariant mass m23.
Double_t p2_Parent_
Momentum of track 2 in parent rest frame.
std::vector< Double_t > mass_
Vector of daughter particles masses.
const Double_t m1Sq_
Mass of particle 1 squared.
Double_t getm3() const
Get the m3 mass.
Double_t getm1() const
Get the m1 mass.
Double_t getm23SqMax() const
Get the m23Sq maximum, (mParent - m1)^2.
Double_t m13_
Invariant mass m13.
Double_t p1_13_
Momentum of track 1 in 1-3 rest frame.
void drawDPContour(Int_t orientation=1323, Int_t nbins=100)
Method to draw the Dalitz plot contours on the top of the histo previously drawn. ...
Bool_t warnings_
Enable/disable warning messages.
Double_t getm13SqMax() const
Get the m13Sq maximum, (mParent - m2)^2.
Double_t distanceFromDPCentre() const
Calculate the distance from the currently set (m13Sq, m23Sq) point to the centre of the Dalitz plot (...
const Double_t mSqDTot_
Sum of the squares of the daughter masses.
Double_t getm12SqMin() const
Get the m12Sq minimum, (m1 + m2)^2.
std::vector< Double_t > mMin_
Vector of the minimum mij values.
File containing declaration of LauKinematics class.
Double_t genm23Sq() const
Randomly generate the invariant mass squared m23Sq.
void calcm12Sq()
Calculate m12Sq from m13Sq and m23Sq.
Double_t p2_23_
Momentum of track 2 in 2-3 rest frame.
Double_t thetaPrime_
theta&#39; co-ordinate
std::vector< Double_t > mSqDiff_
Vector of the difference between the mSqMax and mSqMin.
Double_t calcSqDPJacobian() const
Calculate the Jacobian for the transformation m23^2, m13^2 -&gt; m&#39;, theta&#39; (square DP) at the currently ...
void updateKinematics(Double_t m13Sq, Double_t m23Sq)
Update all kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
Double_t getm2() const
Get the m2 mass.
Double_t p2_13_
Momentum of track 2 in 1-3 rest frame.
void genFlatSqDP(Double_t &mPrime, Double_t &thetaPrime) const
Routine to generate events flat in the square Dalitz plot.
Double_t p1_Parent_
Momentum of track 1 in parent rest frame.
Double_t pCalc(Double_t energy, Double_t massSq) const
General method to calculate the momentum of a particle, given its energy and invariant mass squared...
const Double_t m3Sq_
Mass of particle 3 squared.
Bool_t withinSqDPLimits(Double_t mPrime, Double_t thetaPrime) const
Check whether a given (m&#39;,theta&#39;) point is within the kinematic limits of the Dalitz plot...
Double_t mFromC(Double_t mijSq, Double_t cij, Double_t mij, Int_t i, Int_t j, Int_t k) const
General method to calculate mikSq given mijSq and cosHel_ij.
void flipAndUpdateKinematics()
Flips the DP variables m13^2 &lt;-&gt; m23^2 and recalculates all kinematic quantities. ...
void updateMassSq_m23(Double_t m23, Double_t c23)
Update the variables m12Sq_ and m13Sq_ given the invariant mass m23 and the cosine of the helicity an...
const Double_t invPi
One over Pi.
const Double_t pi
Pi.
Definition: LauConstants.hh:89
Double_t getmParent() const
Get parent mass.
void updateSqDPMassSquares(Double_t mPrime, Double_t thetaPrime)
Update some kinematic quantities based on the square DP co-ordinates m&#39; and theta&#39;.
Double_t m13Sq_
Invariant mass m13 squared.
Bool_t withinDPLimits(Double_t m13Sq, Double_t m23Sq) const
Check whether a given (m13Sq,m23Sq) point is within the kinematic limits of the Dalitz plot...
Double_t m12Sq_
Invariant mass m12 squared.
Double_t qk_
Momentum q of particle k.
File containing LauRandom namespace.
Double_t p3_Parent_
Momentum of track 3 in parent rest frame.
void updateMassSq_m12(Double_t m12, Double_t c12)
Update the variables m23Sq_ and m13Sq_ given the invariant mass m12 and the cosine of the helicity an...
void rotateAndUpdateKinematics()
Cyclically rotates the DP variables (m12 -&gt; m23, m23 -&gt; m13, m13 -&gt; m12) and recalculates all kinemat...
Double_t genm13Sq() const
Randomly generate the invariant mass squared m13Sq.
std::vector< Double_t > mSqMax_
Vector of the maximum mijSq values.
Double_t c23_
Cosine of the helicity angle theta23, which is defined as the angle between 1&amp;2 in the rest frame of ...
Double_t dal(Double_t *x, Double_t *par)
void calcParentFrameMomenta()
Calculate the momenta of each daughter in the parent rest frame.
std::vector< Double_t > mDiff_
Vector of the difference between the mMax and mMin.
Bool_t squareDP_
Should we calculate the square DP co-ordinates or not?
File containing LauConstants namespace.
Double_t m12_
Invariant mass m12.
Double_t c12_
Cosine of the helicity angle theta12, which is defined as the angle between 1&amp;3 in the rest frame of ...
const Double_t mParent_
Mass of parent particle.
Double_t p3_12_
Momentum of track 3 in 1-2 rest frame.
Double_t m23Sq_
Invariant mass m23 squared.
Double_t getm12SqMax() const
Get the m12Sq maximum, (mParent - m3)^2.
void genFlatPhaseSpace(Double_t &m13Sq, Double_t &m23Sq) const
Routine to generate events flat in phase-space.
virtual ~LauKinematics()
Destructor.
Class for calculating 3-body kinematic quantities.
Double_t value() const
The value of the parameter.
Double_t genm12Sq() const
Randomly generate the invariant mass squared m12Sq.
const Double_t piBy2
Pi divided by two.
Double_t c13_
Cosine of the helicity angle theta13, which is defined as the angle between 1&amp;2 in the rest frame of ...
void calcHelicities()
Calculate cosines of the helicity angles, momenta of daughters and bachelor in various ij rest frames...
const Double_t m2Sq_
Mass of particle 2 squared.
Double_t cFromM(Double_t mijSq, Double_t mikSq, Double_t mij, Int_t i, Int_t j, Int_t k) const
General method to calculate the cos(helicity) variables from the masses of the particles.
Double_t p1_12_
Momentum of track 1 in 1-2 rest frame.
Double_t mPrime_
m&#39; co-ordinate
void updateMassSquares(Double_t m13Sq, Double_t m23Sq)
Update some kinematic quantities based on the DP co-ordinates m13Sq and m23Sq.
std::vector< Double_t > mSq_
Vector of daughter particles masses squared.
void calcSqDPVars()
Calculate the m&#39; and theta&#39; variables for the square Dalitz plot.