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