laura is hosted by Hepforge, IPPP Durham
Laura++  3.6.0
A maximum likelihood fitting package for performing Dalitz-plot analysis.
LauDPDepCruijffPdf.cc
Go to the documentation of this file.
1 
2 /*
3 Copyright 2009 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 #include <vector>
31 using std::cerr;
32 using std::cout;
33 using std::endl;
34 
35 #include "LauConstants.hh"
36 #include "LauDPDepCruijffPdf.hh"
37 #include "LauDaughters.hh"
38 #include "LauKinematics.hh"
39 
40 #include "TMath.h"
41 #include "TSystem.h"
42 
43 LauDPDepCruijffPdf::LauDPDepCruijffPdf( const TString& theVarName,
44  const std::vector<LauAbsRValue*>& params,
45  Double_t minAbscissa,
46  Double_t maxAbscissa,
47  const LauDaughters* daughters,
48  const std::vector<Double_t>& meanCoeffs,
49  const std::vector<Double_t>& sigmaLCoeffs,
50  const std::vector<Double_t>& sigmaRCoeffs,
51  const std::vector<Double_t>& alphaLCoeffs,
52  const std::vector<Double_t>& alphaRCoeffs,
53  DPAxis dpAxis ) :
54  LauAbsPdf( theVarName, params, minAbscissa, maxAbscissa ),
55  kinematics_( daughters ? daughters->getKinematics() : 0 ),
56  mean_( 0 ),
57  sigmaL_( 0 ),
58  sigmaR_( 0 ),
59  alphaL_( 0 ),
60  alphaR_( 0 ),
61  meanVal_( 0.0 ),
62  sigmaLVal_( 0.0 ),
63  sigmaRVal_( 0.0 ),
64  alphaLVal_( 0.0 ),
65  alphaRVal_( 0.0 ),
66  meanCoeffs_( meanCoeffs ),
67  sigmaLCoeffs_( sigmaLCoeffs ),
68  sigmaRCoeffs_( sigmaRCoeffs ),
69  alphaLCoeffs_( alphaLCoeffs ),
70  alphaRCoeffs_( alphaRCoeffs ),
71  dpAxis_( dpAxis )
72 {
73  // Constructor for the Dalitz-plot dependent Cruijff PDF.
74 
75  // Check we have a valid kinematics object
76  if ( ! kinematics_ ) {
77  cerr << "ERROR in LauDPDepCruijffPdf::LauDPDepCruijffPdf : Have not been provided with a valid DP kinematics object."
78  << endl;
79  gSystem->Exit( EXIT_FAILURE );
80  }
81 
82  // The parameters in params are the mean, sigmaR, sigmaL, alphaR and alphaL.
83  // The last two arguments specify the range in which the PDF is defined, and the PDF
84  // will be normalised w.r.t. these limits.
85 
86  mean_ = this->findParameter( "mean" );
87  sigmaR_ = this->findParameter( "sigmaR" );
88  sigmaL_ = this->findParameter( "sigmaL" );
89  alphaR_ = this->findParameter( "alphaR" );
90  alphaL_ = this->findParameter( "alphaL" );
91 
92  if ( ( this->nParameters() != 5 ) || ( mean_ == 0 ) || ( sigmaR_ == 0 ) || ( sigmaL_ == 0 ) ||
93  ( alphaL_ == 0 ) || ( alphaR_ == 0 ) ) {
94  cerr << "ERROR in LauDPDepCruijffPdf constructor: LauDPDepCruijffPdf requires 5 parameters: \"mean\", \"sigmaL\", \"sigmaR\", \"alphaR\" and \"alphaL\"."
95  << endl;
96  gSystem->Exit( EXIT_FAILURE );
97  }
98 
99  // Cache the normalisation factor
100  this->calcNorm();
101 }
102 
104 {
105  // Destructor
106 }
107 
109 {
110  this->setNorm( 1.0 );
111 }
112 
114 {
115  // Check that the given abscissa is within the allowed range
116  if ( ! this->checkRange( abscissas ) ) {
117  gSystem->Exit( EXIT_FAILURE );
118  }
119 
120  // Get our abscissa
121  Double_t abscissa = abscissas[0];
122 
123  // Get the up to date parameter values
129 
130  // Find out the DP position
131  Double_t dpPos( 0.0 );
132  UInt_t nVars = this->nInputVars();
133  if ( abscissas.size() == nVars + 2 ) {
134  Double_t m13Sq = abscissas[nVars];
135  Double_t m23Sq = abscissas[nVars + 1];
136 
137  if ( dpAxis_ == M12 ) {
138  dpPos = kinematics_->calcThirdMassSq( m13Sq, m23Sq );
139  } else if ( dpAxis_ == M13 ) {
140  dpPos = m13Sq;
141  } else if ( dpAxis_ == M23 ) {
142  dpPos = m23Sq;
143  } else if ( dpAxis_ == MMIN ) {
144  dpPos = TMath::Min( m13Sq, m23Sq );
145  } else if ( dpAxis_ == MMAX ) {
146  dpPos = TMath::Max( m13Sq, m23Sq );
147  } else {
148  dpPos = kinematics_->distanceFromDPCentre( m13Sq, m23Sq );
149  }
150  }
151 
152  // Scale the parameters according to the DP position
153  this->scalePars( dpPos );
154 
155  // Calculate the PDF value
156  Double_t value = this->currentPDFValue( abscissa );
157 
158  // Calculate the normalisation
159  IntMethod sumMethod = this->integMethod();
160  Double_t normFac = ( sumMethod == GaussLegendre ) ? this->integrGaussLegendre()
161  : this->integTrapezoid();
162 
163  value /= normFac;
164 
165  this->setUnNormPDFVal( value );
166 }
167 
168 Double_t LauDPDepCruijffPdf::currentPDFValue( Double_t abscissa ) const
169 {
170  Double_t arg = abscissa - meanVal_;
171  Double_t coef( 0.0 );
172  if ( arg < 0.0 ) {
173  if ( TMath::Abs( sigmaLVal_ ) > 1e-30 ) {
174  coef = -1.0 / ( 2.0 * sigmaLVal_ * sigmaLVal_ + alphaLVal_ * arg * arg );
175  }
176  } else {
177  if ( TMath::Abs( sigmaRVal_ ) > 1e-30 ) {
178  coef = -1.0 / ( 2.0 * sigmaRVal_ * sigmaRVal_ + alphaRVal_ * arg * arg );
179  }
180  }
181  return TMath::Exp( coef * arg * arg );
182 }
183 
185 {
186  if ( ! this->normWeightsDone() ) {
187  this->getNormWeights();
188  }
189 
190  const std::vector<LauAbscissas>& norm_abscissas = this->normAbscissas();
191  const std::vector<Double_t>& norm_weights = this->normWeights();
192 
193  // Now compute the integral
194  Double_t norm( 0.0 );
195  for ( UInt_t i = 0; i < norm_weights.size(); i++ ) {
196  Double_t fun = this->currentPDFValue( norm_abscissas[i][0] );
197  Double_t intFactor = 0.5 * this->getRange();
198  norm += norm_weights[i] * intFactor * fun;
199  }
200 
201  return norm;
202 }
203 
205 {
206  static Double_t norm( 0.0 );
207  Double_t range = this->getRange();
208 
209  if ( this->nNormPoints() == 1 ) {
210 
211  Double_t abscissa = this->getMinAbscissa();
212  Double_t funAbsMin = this->currentPDFValue( abscissa );
213 
214  abscissa = this->getMaxAbscissa();
215  Double_t funAbsMax = this->currentPDFValue( abscissa );
216 
217  norm = 0.5 * range * ( funAbsMin + funAbsMax );
218  return norm;
219 
220  } else {
221 
222  Double_t abscVal( 0.0 ), tnm( 0.0 ), sum( 0.0 ), del( 0.0 );
223  Int_t it( 0 ), j( 0 );
224 
225  for ( it = 1, j = 1; j < this->nNormPoints() - 1; j++ ) {
226  it <<= 1;
227  }
228  tnm = it;
229  del = range / tnm;
230  abscVal = this->getMinAbscissa() + 0.5 * del;
231 
232  for ( sum = 0.0, j = 1; j < it; j++, abscVal += del ) {
233 
234  Double_t funVal = this->currentPDFValue( abscVal );
235  sum += funVal;
236  }
237 
238  norm = 0.5 * ( norm + sum * range / tnm );
239  return norm;
240  }
241 }
242 
243 void LauDPDepCruijffPdf::scalePars( Double_t dpPos )
244 {
245  Int_t power = 1;
246  for ( std::vector<Double_t>::const_iterator iter = meanCoeffs_.begin(); iter != meanCoeffs_.end();
247  ++iter ) {
248  Double_t coeff = ( *iter );
249  meanVal_ += coeff * TMath::Power( dpPos, power );
250  ++power;
251  }
252 
253  power = 1;
254  for ( std::vector<Double_t>::const_iterator iter = sigmaLCoeffs_.begin();
255  iter != sigmaLCoeffs_.end();
256  ++iter ) {
257  Double_t coeff = ( *iter );
258  sigmaLVal_ += coeff * TMath::Power( dpPos, power );
259  ++power;
260  }
261 
262  power = 1;
263  for ( std::vector<Double_t>::const_iterator iter = sigmaRCoeffs_.begin();
264  iter != sigmaRCoeffs_.end();
265  ++iter ) {
266  Double_t coeff = ( *iter );
267  sigmaRVal_ += coeff * TMath::Power( dpPos, power );
268  ++power;
269  }
270 
271  power = 1;
272  for ( std::vector<Double_t>::const_iterator iter = alphaLCoeffs_.begin();
273  iter != alphaLCoeffs_.end();
274  ++iter ) {
275  Double_t coeff = ( *iter );
276  alphaLVal_ += coeff * TMath::Power( dpPos, power );
277  ++power;
278  }
279 
280  power = 1;
281  for ( std::vector<Double_t>::const_iterator iter = alphaRCoeffs_.begin();
282  iter != alphaRCoeffs_.end();
283  ++iter ) {
284  Double_t coeff = ( *iter );
285  alphaRVal_ += coeff * TMath::Power( dpPos, power );
286  ++power;
287  }
288 }
289 
291 {
292  // Get the up to date parameter values
298 
299  // Find out the DP position
300  Double_t dpPos( 0.0 );
301  if ( dpAxis_ == M12 ) {
302  dpPos = kinematics->getm12Sq();
303  } else if ( dpAxis_ == M13 ) {
304  dpPos = kinematics->getm13Sq();
305  } else if ( dpAxis_ == M23 ) {
306  dpPos = kinematics->getm23Sq();
307  } else {
308  dpPos = kinematics->distanceFromDPCentre();
309  }
310 
311  // Scale the parameters according to the DP position
312  this->scalePars( dpPos );
313 
314  // Calculate the PDF height
315 
316  LauAbscissas maxPoint( 3 );
317  maxPoint[0] = meanVal_;
318  maxPoint[1] = kinematics->getm13Sq();
319  maxPoint[2] = kinematics->getm23Sq();
320 
321  if ( meanVal_ < this->getMinAbscissa() ) {
322  maxPoint[0] = this->getMinAbscissa();
323  } else if ( meanVal_ > this->getMaxAbscissa() ) {
324  maxPoint[0] = this->getMaxAbscissa();
325  }
326 
327  this->calcLikelihoodInfo( maxPoint );
328  Double_t height = this->getUnNormLikelihood();
329 
330  // Multiply by a small factor to avoid problems from rounding errors
331  height *= 1.01;
332 
333  this->setMaxHeight( height );
334 }
DPAxis
Define possibilties for the DP axes.
LauAbsRValue * alphaL_
Alpha of left Gaussian.
Double_t sigmaRVal_
Sigma of right Gaussian.
const std::vector< Double_t > sigmaRCoeffs_
Coefficients of sigma for the right Gaussian.
virtual Double_t getMaxAbscissa() const
Retrieve the maximum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:141
void scalePars(Double_t dpPos)
Scale parameters by their dependence on the DP position.
LauAbsRValue * alphaR_
Alpha of right Gaussian.
virtual void calcLikelihoodInfo(const LauAbscissas &abscissas)
Calculate the likelihood (and intermediate info) for a given abscissa.
Double_t value() const
The value of the parameter.
virtual Double_t getMinAbscissa() const
Retrieve the minimum value of the (primary) abscissa.
Definition: LauAbsPdf.hh:135
Double_t getm12Sq() const
Get the m12 invariant mass square.
Double_t distanceFromDPCentre() const
Calculate the distance from the currently set (m13Sq, m23Sq) point to the centre of the Dalitz plot (...
LauAbsRValue * sigmaL_
Sigma of left Gaussian.
IntMethod
The possible numerical intergration methods.
Definition: LauAbsPdf.hh:61
std::vector< Double_t > LauAbscissas
The type used for containing multiple abscissa values.
Definition: LauAbsPdf.hh:58
const std::vector< Double_t > sigmaLCoeffs_
Coefficients of sigma for the left Gaussian.
virtual Double_t getRange() const
Retrieve the range of the (primary) abscissa.
Definition: LauAbsPdf.hh:147
Double_t currentPDFValue(Double_t abscissa) const
Current PDF value.
virtual LauAbsRValue * findParameter(const TString &parName)
Retrieve the specified parameter.
Definition: LauAbsPdf.cc:431
virtual Bool_t normWeightsDone() const
Check whether the normalisation weights have been calculated.
Definition: LauAbsPdf.hh:472
Double_t meanVal_
Gaussian mean.
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)
virtual void calcNorm()
Calculate the normalisation.
LauAbsRValue * mean_
Gaussian mean.
const LauKinematics * kinematics_
The current DP kinematics.
virtual UInt_t nParameters() const
Retrieve the number of PDF parameters.
Definition: LauAbsPdf.hh:109
virtual void setNorm(Double_t norm)
Set the normalisation factor.
Definition: LauAbsPdf.hh:348
virtual Double_t unblindValue() const =0
The unblinded value of the parameter.
Double_t alphaRVal_
Alpha of right Gaussian.
Double_t alphaLVal_
Alpha of left Gaussian.
virtual UInt_t nInputVars() const
Retrieve the number of abscissas.
Definition: LauAbsPdf.hh:121
File containing declaration of LauDaughters class.
LauAbsRValue * sigmaR_
Sigma of right Gaussian.
virtual void setUnNormPDFVal(Double_t unNormPDFVal)
Set the unnormalised likelihood.
Definition: LauAbsPdf.hh:394
virtual ~LauDPDepCruijffPdf()
Destructor.
File containing declaration of LauDPDepCruijffPdf class.
Double_t sigmaLVal_
Sigma of left Gaussian.
Class for defining the abstract interface for PDF classes.
Definition: LauAbsPdf.hh:54
LauDPDepCruijffPdf(const TString &theVarName, const std::vector< LauAbsRValue * > &params, Double_t minAbscissa, Double_t maxAbscissa, const LauDaughters *daughters, const std::vector< Double_t > &meanCoeffs, const std::vector< Double_t > &sigmaLCoeffs, const std::vector< Double_t > &sigmaRCoeffs, const std::vector< Double_t > &alphaLCoeffs, const std::vector< Double_t > &alphaRCoeffs, DPAxis dpAxis)
Constructor.
virtual const std::vector< Double_t > & normWeights() const
Retrieve the weights used for normalisation.
Definition: LauAbsPdf.hh:493
const std::vector< Double_t > alphaRCoeffs_
Coefficients of alpha for the right Gaussian.
File containing LauConstants namespace.
DPAxis dpAxis_
The DP axis we depend on.
virtual void calcPDFHeight(const LauKinematics *kinematics)
Calculate the PDF height.
virtual void getNormWeights()
Calculate the weights and abscissas used for normalisation.
Definition: LauAbsPdf.cc:527
virtual void setMaxHeight(Double_t maxHeight)
Set the maximum height.
Definition: LauAbsPdf.hh:354
virtual const std::vector< LauAbscissas > & normAbscissas() const
Retrieve the abscissa points used for normalisation.
Definition: LauAbsPdf.hh:487
Double_t range() const
The range allowed for the parameter.
const std::vector< Double_t > meanCoeffs_
Coefficients of Gaussian mean.
Class for calculating 3-body kinematic quantities.
virtual IntMethod integMethod() const
Retrieve the integration method used to normalise the PDF.
Definition: LauAbsPdf.hh:317
virtual Double_t integTrapezoid()
Integrate the PDF using the simple trapezoid method.
Class that defines the particular 3-body decay under study.
Definition: LauDaughters.hh:47
Double_t getm23Sq() const
Get the m23 invariant mass square.
virtual Int_t nNormPoints() const
Retrieve the number of points to integrate over when normalising.
Definition: LauAbsPdf.hh:305
virtual Double_t getUnNormLikelihood() const
Retrieve the unnormalised likelihood value.
Definition: LauAbsPdf.hh:218
Double_t getm13Sq() const
Get the m13 invariant mass square.
@ GaussLegendre
Definition: LauAbsPdf.hh:62
virtual Double_t integrGaussLegendre()
Integrate the PDF using the Gauss-Legendre method.
File containing declaration of LauKinematics class.
virtual Bool_t checkRange(const LauAbscissas &abscissas) const
Check that all abscissas are within their allowed ranges.
Definition: LauAbsPdf.cc:243
const std::vector< Double_t > alphaLCoeffs_
Coefficients of alpha for the left Gaussian.