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