Logo ROOT   6.08/07
Reference Guide
RooGExpModel.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 /**
18 \file RooGExpModel.cxx
19 \class RooGExpModel
20 \ingroup Roofit
21 
22 Class RooGExpModel is a RooResolutionModel implementation that models
23 a resolution function that is the convolution of a Gaussian with
24 a one-sided exponential. Object of class RooGExpModel can be used
25 for analytical convolutions with classes inheriting from RooAbsAnaConvPdf
26 **/
27 
28 #include "RooFit.h"
29 
30 #include "Riostream.h"
31 #include "Riostream.h"
32 #include "RooGExpModel.h"
33 #include "RooMath.h"
34 #include "RooRealConstant.h"
35 #include "RooRandom.h"
36 #include "RooMath.h"
37 #include "TMath.h"
38 
39 #include "TError.h"
40 
41 using namespace std;
42 
44 ;
45 
46 
47 
48 ////////////////////////////////////////////////////////////////////////////////
49 
50 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
51  RooAbsReal& _sigma, RooAbsReal& _rlife,
52  Bool_t nlo, Type type) :
53  RooResolutionModel(name,title,xIn),
54  sigma("sigma","Width",this,_sigma),
55  rlife("rlife","Life time",this,_rlife),
56  ssf("ssf","Sigma Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
57  rsf("rsf","RLife Scale Factor",this,(RooRealVar&)RooRealConstant::value(1)),
58  _flip(type==Flipped),_nlo(nlo), _flatSFInt(kFALSE), _asympInt(kFALSE)
59 {
60 }
61 
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 
66 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
67  RooAbsReal& _sigma, RooAbsReal& _rlife,
68  RooAbsReal& _rsSF,
69  Bool_t nlo, Type type) :
70  RooResolutionModel(name,title,xIn),
71  sigma("sigma","Width",this,_sigma),
72  rlife("rlife","Life time",this,_rlife),
73  ssf("ssf","Sigma Scale Factor",this,_rsSF),
74  rsf("rsf","RLife Scale Factor",this,_rsSF),
75  _flip(type==Flipped),
76  _nlo(nlo),
79 {
80 }
81 
82 
83 
84 ////////////////////////////////////////////////////////////////////////////////
85 
86 RooGExpModel::RooGExpModel(const char *name, const char *title, RooRealVar& xIn,
87  RooAbsReal& _sigma, RooAbsReal& _rlife,
88  RooAbsReal& _sigmaSF, RooAbsReal& _rlifeSF,
89  Bool_t nlo, Type type) :
90  RooResolutionModel(name,title,xIn),
91  sigma("sigma","Width",this,_sigma),
92  rlife("rlife","Life time",this,_rlife),
93  ssf("ssf","Sigma Scale Factor",this,_sigmaSF),
94  rsf("rsf","RLife Scale Factor",this,_rlifeSF),
95  _flip(type==Flipped),
96  _nlo(nlo),
99 {
100 }
101 
102 
103 
104 ////////////////////////////////////////////////////////////////////////////////
105 
106 RooGExpModel::RooGExpModel(const RooGExpModel& other, const char* name) :
107  RooResolutionModel(other,name),
108  sigma("sigma",this,other.sigma),
109  rlife("rlife",this,other.rlife),
110  ssf("ssf",this,other.ssf),
111  rsf("rsf",this,other.rsf),
112  _flip(other._flip),
113  _nlo(other._nlo),
114  _flatSFInt(other._flatSFInt),
115  _asympInt(other._asympInt)
116 {
117 }
118 
119 
120 
121 ////////////////////////////////////////////////////////////////////////////////
122 /// Destructor
123 
125 {
126 }
127 
128 
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 
132 Int_t RooGExpModel::basisCode(const char* name) const
133 {
134  if (!TString("exp(-@0/@1)").CompareTo(name)) return expBasisPlus ;
135  if (!TString("exp(@0/@1)").CompareTo(name)) return expBasisMinus ;
136  if (!TString("exp(-abs(@0)/@1)").CompareTo(name)) return expBasisSum ;
137  if (!TString("exp(-@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisPlus ;
138  if (!TString("exp(@0/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisMinus ;
139  if (!TString("exp(-abs(@0)/@1)*sin(@0*@2)").CompareTo(name)) return sinBasisSum ;
140  if (!TString("exp(-@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisPlus ;
141  if (!TString("exp(@0/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisMinus ;
142  if (!TString("exp(-abs(@0)/@1)*cos(@0*@2)").CompareTo(name)) return cosBasisSum ;
143  if (!TString("exp(-@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisPlus;
144  if (!TString("exp(@0/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisMinus;
145  if (!TString("exp(-abs(@0)/@1)*sinh(@0*@2/2)").CompareTo(name)) return sinhBasisSum;
146  if (!TString("exp(-@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisPlus;
147  if (!TString("exp(@0/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisMinus;
148  if (!TString("exp(-abs(@0)/@1)*cosh(@0*@2/2)").CompareTo(name)) return coshBasisSum;
149  return 0 ;
150 }
151 
152 
153 
154 ////////////////////////////////////////////////////////////////////////////////
155 
157 {
158  static Double_t root2(sqrt(2.)) ;
159 // static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
160 // static Double_t rootpi(sqrt(atan2(0.,-1.)));
161 
162  BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
163  BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
164 
165  Double_t fsign = _flip?-1:1 ;
166 
167  Double_t sig = sigma*ssf ;
168  Double_t rtau = rlife*rsf ;
169 
171  // added, FMV 07/27/03
172  if (basisType == coshBasis && _basisCode!=noBasis ) {
173  Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
174  if (dGamma==0) basisType = expBasis;
175  }
176 
177  // *** 1st form: Straight GExp, used for unconvoluted PDF or expBasis with 0 lifetime ***
178  if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
179  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 1st form" << endl ;
180 
181  Double_t expArg = sig*sig/(2*rtau*rtau) + fsign*x/rtau ;
182 
183  Double_t result ;
184  if (expArg<300) {
185  result = 1/(2*rtau) * exp(expArg) * RooMath::erfc(sig/(root2*rtau) + fsign*x/(root2*sig));
186  } else {
187  // If exponent argument is very large, bring canceling RooMath::erfc() term inside exponent
188  // to avoid floating point over/underflows of intermediate calculations
189  result = 1/(2*rtau) * exp(expArg + logErfC(sig/(root2*rtau) + fsign*x/(root2*sig))) ;
190  }
191 
192 // Double_t result = 1/(2*rtau)
193 // * exp(sig*sig/(2*rtau*rtau) + fsign*x/rtau)
194 // * RooMath::erfc(sig/(root2*rtau) + fsign*x/(root2*sig));
195 
196  // equivalent form, added FMV, 07/24/03
197  //Double_t xprime = x/rtau ;
198  //Double_t c = sig/(root2*rtau) ;
199  //Double_t u = xprime/(2*c) ;
200  //Double_t result = 0.5*evalCerf(fsign*u,c).real() ; // sign=-1 !
201 
202  if (_basisCode!=0 && basisSign==Both) result *= 2 ;
203  //cout << "1st form " << "x= " << x << " result= " << result << endl;
204  return result ;
205  }
206 
207  // *** 2nd form: 0, used for sinBasis and cosBasis with tau=0 ***
208  if (tau==0) {
209  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 2nd form" << endl ;
210  return 0. ;
211  }
212 
213  Double_t omega = (basisType!=expBasis)?((RooAbsReal*)basis().getParameter(2))->getVal():0. ;
214 
215  // *** 3nd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
216  if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
217  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 3d form tau=" << tau << endl ;
218  Double_t result(0) ;
219  if (basisSign!=Minus) result += calcDecayConv(+1,tau,sig,rtau,fsign) ; // modified FMV,08/13/03
220  if (basisSign!=Plus) result += calcDecayConv(-1,tau,sig,rtau,fsign) ; // modified FMV,08/13/03
221  //cout << "3rd form " << "x= " << x << " result= " << result << endl;
222  return result ;
223  }
224 
225  // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
226  Double_t wt = omega *tau ;
227  if (basisType==sinBasis) {
228  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName() << ") 4th form omega = "
229  << omega << ", tau = " << tau << endl ;
230  Double_t result(0) ;
231  if (wt==0.) return result ;
232  if (basisSign!=Minus) result += -1*calcSinConv(+1,sig,tau,omega,rtau,fsign).imag() ;
233  if (basisSign!=Plus) result += -1*calcSinConv(-1,sig,tau,omega,rtau,fsign).imag() ;
234  //cout << "4th form " << "x= " << x << " result= " << result << endl;
235  return result ;
236  }
237 
238  // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
239  if (basisType==cosBasis) {
240  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
241  << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
242  Double_t result(0) ;
243  if (basisSign!=Minus) result += calcSinConv(+1,sig,tau,omega,rtau,fsign).real() ;
244  if (basisSign!=Plus) result += calcSinConv(-1,sig,tau,omega,rtau,fsign).real() ;
245  //cout << "5th form " << "x= " << x << " result= " << result << endl;
246  return result ;
247  }
248 
249 
250  // *** 6th form: Convolution with exp(-t/tau)*sinh(dgamma*t/2), used for sinhBasis ***
251  if (basisType==sinhBasis) {
252  Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
253 
254  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
255  << ") 6th form = " << dgamma << ", tau = " << tau << endl;
256  Double_t result(0);
257  //if (basisSign!=Minus) result += calcSinhConv(+1,+1,-1,tau,dgamma,sig,rtau,fsign);
258  //if (basisSign!=Plus) result += calcSinhConv(-1,-1,+1,tau,dgamma,sig,rtau,fsign);
259  // better form, since it also accounts for the numerical divergence region, added FMV, 07/24/03
260  Double_t tau1 = 1/(1/tau-dgamma/2) ;
261  Double_t tau2 = 1/(1/tau+dgamma/2) ;
262  if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)-calcDecayConv(+1,tau2,sig,rtau,fsign));
263  // modified FMV,08/13/03
264  if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau2,sig,rtau,fsign)-calcDecayConv(-1,tau1,sig,rtau,fsign));
265  // modified FMV,08/13/03
266  //cout << "6th form " << "x= " << x << " result= " << result << endl;
267  return result;
268  }
269 
270  // *** 7th form: Convolution with exp(-t/tau)*cosh(dgamma*t/2), used for coshBasis ***
271  if (basisType==coshBasis) {
272  Double_t dgamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
273 
274  if (verboseEval()>2) cout << "RooGExpModel::evaluate(" << GetName()
275  << ") 7th form = " << dgamma << ", tau = " << tau << endl;
276  Double_t result(0);
277  //if (basisSign!=Minus) result += calcCoshConv(+1,tau,dgamma,sig,rtau,fsign);
278  //if (basisSign!=Plus) result += calcCoshConv(-1,tau,dgamma,sig,rtau,fsign);
279  // better form, since it also accounts for the numerical divergence region, added FMV, 07/24/03
280  Double_t tau1 = 1/(1/tau-dgamma/2) ;
281  Double_t tau2 = 1/(1/tau+dgamma/2) ;
282  if (basisSign!=Minus) result += 0.5*(calcDecayConv(+1,tau1,sig,rtau,fsign)+calcDecayConv(+1,tau2,sig,rtau,fsign));
283  // modified FMV,08/13/03
284  if (basisSign!=Plus) result += 0.5*(calcDecayConv(-1,tau1,sig,rtau,fsign)+calcDecayConv(-1,tau2,sig,rtau,fsign));
285  // modified FMV,08/13/03
286  //cout << "7th form " << "x= " << x << " result= " << result << endl;
287  return result;
288  }
289  R__ASSERT(0) ;
290  return 0 ;
291  }
292 
293 
294 
295 ////////////////////////////////////////////////////////////////////////////////
296 /// Approximation of the log of the complex error function
297 
299 {
300  Double_t t,z,ans;
301  z=fabs(xx);
302  t=1.0/(1.0+0.5*z);
303 
304  if(xx >= 0.0)
305  ans=log(t)+(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
306  t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277)))))))));
307  else
308  ans=log(2.0-t*exp(-z*z-1.26551223+t*(1.00002368+t*(0.37409196+t*(0.09678418+t*(-0.18628806+
309  t*(0.27886807+t*(-1.13520398+t*(1.48851587+t*(-0.82215223+t*0.17087277))))))))));
310 
311  return ans;
312 }
313 
314 
315 
316 ////////////////////////////////////////////////////////////////////////////////
317 
318 std::complex<Double_t> RooGExpModel::calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t omega, Double_t rtau, Double_t fsign) const
319 {
320  static Double_t root2(sqrt(2.)) ;
321 
322  Double_t s1= -sign*x/tau;
323  //Double_t s1= x/tau;
324  Double_t c1= sig/(root2*tau);
325  Double_t u1= s1/(2*c1);
326  Double_t s2= x/rtau;
327  Double_t c2= sig/(root2*rtau);
328  Double_t u2= fsign*s2/(2*c2) ;
329  //Double_t u2= s2/(2*c2) ;
330 
331  std::complex<Double_t> eins(1,0);
332  std::complex<Double_t> k(1/tau,sign*omega);
333  //return (evalCerf(-sign*omega*tau,u1,c1)+evalCerf(0,u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
334 
335  return (evalCerf(-sign*omega*tau,u1,c1)+std::complex<Double_t>(evalCerfRe(u2,c2),0)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
336  // equivalent form, added FMV, 07/24/03
337  //return (evalCerf(-sign*omega*tau,-sign*u1,c1)+evalCerf(0,fsign*u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
338 }
339 
340 
341 // added FMV,08/18/03
342 
343 ////////////////////////////////////////////////////////////////////////////////
344 
346 {
347  static Double_t root2(sqrt(2.)) ;
348 
349  Double_t s1= -sign*x/tau;
350  //Double_t s1= x/tau;
351  Double_t c1= sig/(root2*tau);
352  Double_t u1= s1/(2*c1);
353  Double_t s2= x/rtau;
354  Double_t c2= sig/(root2*rtau);
355  Double_t u2= fsign*s2/(2*c2) ;
356  //Double_t u2= s2/(2*c2) ;
357 
358  Double_t eins(1);
359  Double_t k(1/tau);
360  return (evalCerfRe(u1,c1)+evalCerfRe(u2,c2)*fsign*sign) / (eins + k*fsign*sign*rtau) ;
361  // equivalent form, added FMV, 07/24/03
362  //return (evalCerf(-sign*u1,c1).real()+evalCerf(fsign*u2,c2).real()*fsign*sign) / (eins + k*fsign*sign*rtau) ;
363 }
364 
365 
366 
367 ////////////////////////////////////////////////////////////////////////////////
368 
370 // modified FMV,08/13/03
371 {
372  static Double_t root2(sqrt(2.)) ;
373  static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
374  static Double_t rootpi(sqrt(atan2(0.,-1.)));
375 
376  // Process flip status
377  Double_t xp(x) ;
378  //if (_flip) {
379  // xp *= -1 ;
380  // sign *= -1 ;
381  //}
382  xp *= fsign ; // modified FMV,08/13/03
383  sign *= fsign ; // modified FMV,08/13/03
384 
385  Double_t cFly;
386  if ((sign<0)&&(fabs(tau-rtau)<tau/260)) {
387 
388  Double_t MeanTau=0.5*(tau+rtau);
389  if (fabs(xp/MeanTau)>300) {
390  return 0 ;
391  }
392 
393  cFly=1./(MeanTau*MeanTau*root2pi) *
394  exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
395  *(sig*exp(-1/(2*sig*sig)*TMath::Power((sig*sig/MeanTau+xp),2))
396  -(sig*sig/MeanTau+xp)*(rootpi/root2)*RooMath::erfc(sig/(root2*MeanTau)+xp/(root2*sig)));
397 
398  if(_nlo) {
399  Double_t epsilon=0.5*(tau-rtau);
400  Double_t a=sig/(root2*MeanTau)+xp/(root2*sig);
401  cFly += 1./(MeanTau*MeanTau)
402  *exp(-(-xp/MeanTau-sig*sig/(2*MeanTau*MeanTau)))
403  *0.5/MeanTau*epsilon*epsilon*
404  (exp(-a*a)*(sig/MeanTau*root2/rootpi
405  -(4*a*sig*sig)/(2*rootpi*MeanTau*MeanTau)
406  +(-4/rootpi+8*a*a/rootpi)/6
407  *TMath::Power(sig/(root2*MeanTau),3)
408  +2/rootpi*(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
409  (sig/(root2*MeanTau)-a*(sig*sig)/(2*MeanTau*MeanTau))
410  +2/rootpi*((3*sig*sig)/(2*MeanTau*MeanTau)+xp/MeanTau+
411  0.5*TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,2))*sig/(root2*MeanTau))
412  -(2*sig*sig/(MeanTau*MeanTau)+xp/MeanTau+(sig*sig/(MeanTau*MeanTau)+xp/MeanTau)*
413  (3*sig*sig/(2*MeanTau*MeanTau)+xp/MeanTau)
414  +TMath::Power(sig*sig/(MeanTau*MeanTau)+xp/MeanTau,3)/6)*RooMath::erfc(a));
415  }
416 
417  } else {
418 
419  Double_t expArg1 = sig*sig/(2*tau*tau)-sign*xp/tau ;
420  Double_t expArg2 = sig*sig/(2*rtau*rtau)+xp/rtau ;
421 
422  Double_t term1, term2 ;
423  if (expArg1<300) {
424  term1 = exp(expArg1) *RooMath::erfc(sig/(root2*tau)-sign*xp/(root2*sig)) ;
425  } else {
426  term1 = exp(expArg1+logErfC(sig/(root2*tau)-sign*xp/(root2*sig))) ; ;
427  }
428  if (expArg2<300) {
429  term2 = exp(expArg2) *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)) ;
430  } else {
431  term2 = exp(expArg2+logErfC(sig/(root2*rtau)+xp/(root2*sig))) ; ;
432  }
433 
434  cFly=(term1+sign*term2)/(2*(tau+sign*rtau));
435 
436  // WVE prevent numeric underflows
437  if (cFly<1e-100) {
438  cFly = 0 ;
439  }
440 
441  // equivalent form, added FMV, 07/24/03
442  //cFly = calcSinConv(sign, sig, tau, rtau, fsign)/(2*tau) ;
443  }
444 
445  return cFly*2*tau ;
446 }
447 
448 /* commented FMV, 07/24/03
449 
450 ////////////////////////////////////////////////////////////////////////////////
451 
452 Double_t RooGExpModel::calcCoshConv(Double_t sign, Double_t tau, Double_t dgamma, Double_t sig, Double_t rtau, Double_t fsign) const
453 {
454 
455 
456  static Double_t root2(sqrt(2.)) ;
457  static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
458  static Double_t rootpi(sqrt(atan2(0.,-1.)));
459  Double_t tau1 = 1/(1/tau-dgamma/2);
460  Double_t tau2 = 1/(1/tau+dgamma/2);
461  Double_t cFly;
462  Double_t xp(x);
463 
464  //if (_flip) {
465  // xp *= -1 ;
466  // sign *= -1 ;
467  //}
468  xp *= fsign ; // modified FMV,08/13/03
469  sign *= fsign ; // modified FMV,08/13/03
470 
471  cFly=tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
472  *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
473  +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
474  *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
475  +tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
476  *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
477  +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
478  *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));;
479  return cFly;
480 }
481 */
482 
483 /* commented FMV, 07/24/03
484 
485 ////////////////////////////////////////////////////////////////////////////////
486 
487 Double_t RooGExpModel::calcSinhConv(Double_t sign, Double_t sign1, Double_t sign2, Double_t tau, Double_t dgamma, Double_t sig, Double_t rtau, Double_t fsign) const
488 {
489  static Double_t root2(sqrt(2.)) ;
490  static Double_t root2pi(sqrt(2*atan2(0.,-1.))) ;
491  static Double_t rootpi(sqrt(atan2(0.,-1.)));
492  Double_t tau1 = 1/(1/tau-dgamma/2);
493  Double_t tau2 = 1/(1/tau+dgamma/2);
494  Double_t cFly;
495  Double_t xp(x);
496 
497  //if (_flip) {
498  // xp *= -1 ;
499  // sign1 *= -1 ;
500  // sign2 *= -1 ;
501  //}
502  xp *= fsign ; // modified FMV,08/13/03
503  sign1 *= fsign ; // modified FMV,08/13/03
504  sign2 *= fsign ; // modified FMV,08/13/03
505 
506  cFly=sign1*tau1*(exp(sig*sig/(2*tau1*tau1)-sign*xp/tau1)
507  *RooMath::erfc(sig/(root2*tau1)-sign*xp/(root2*sig))
508  +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
509  *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau1+sign*rtau))
510  +sign2*tau2*(exp(sig*sig/(2*tau2*tau2)-sign*xp/tau2)
511  *RooMath::erfc(sig/(root2*tau2)-sign*xp/(root2*sig))
512  +sign*exp(sig*sig/(2*rtau*rtau)+xp/rtau)
513  *RooMath::erfc(sig/(root2*rtau)+xp/(root2*sig)))/(2*(tau2+sign*rtau));;
514  return cFly;
515 }
516 */
517 
518 
519 ////////////////////////////////////////////////////////////////////////////////
520 
521 Int_t RooGExpModel::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
522 {
523  switch(_basisCode) {
524 
525  // Analytical integration capability of raw PDF
526  case noBasis:
527  if (matchArgs(allVars,analVars,convVar())) return 1 ;
528  break ;
529 
530  // Analytical integration capability of convoluted PDF
531  case expBasisPlus:
532  case expBasisMinus:
533  case expBasisSum:
534  case sinBasisPlus:
535  case sinBasisMinus:
536  case sinBasisSum:
537  case cosBasisPlus:
538  case cosBasisMinus:
539  case cosBasisSum:
540  case sinhBasisPlus:
541  case sinhBasisMinus:
542  case sinhBasisSum:
543  case coshBasisPlus:
544  case coshBasisMinus:
545  case coshBasisSum:
546 
547  // Optionally advertise flat integral over sigma scale factor
548  if (_flatSFInt) {
549  if (matchArgs(allVars,analVars,RooArgSet(convVar(),ssf.arg()))) {
550  return 2 ;
551  }
552  }
553 
554  if (matchArgs(allVars,analVars,convVar())) return 1 ;
555  break ;
556  }
557 
558  return 0 ;
559 }
560 
561 
562 
563 ////////////////////////////////////////////////////////////////////////////////
564 
565 Double_t RooGExpModel::analyticalIntegral(Int_t code, const char* rangeName) const
566 {
567  static Double_t root2 = sqrt(2.) ;
568 // static Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
569  Double_t ssfInt(1.0) ;
570 
571  // Code must be 1 or 2
572  R__ASSERT(code==1||code==2) ;
573  if (code==2) {
574  ssfInt = (ssf.max(rangeName)-ssf.min(rangeName)) ;
575  }
576 
577  BasisType basisType = (BasisType)( (_basisCode == 0) ? 0 : (_basisCode/10) + 1 );
578  BasisSign basisSign = (BasisSign)( _basisCode - 10*(basisType-1) - 2 ) ;
579 
581 
582  // added FMV, 07/24/03
583  if (basisType == coshBasis && _basisCode!=noBasis ) {
584  Double_t dGamma = ((RooAbsReal*)basis().getParameter(2))->getVal();
585  if (dGamma==0) basisType = expBasis;
586  }
587  Double_t fsign = _flip?-1:1 ;
588  Double_t sig = sigma*ssf ;
589  Double_t rtau = rlife*rsf ;
590 
591  // *** 1st form????
592  if (basisType==none || ((basisType==expBasis || basisType==cosBasis) && tau==0.)) {
593  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 1st form" << endl ;
594 
595  //Double_t result = 1.0 ; // WVE inferred from limit(tau->0) of cosBasisNorm
596  // finite+asymtotic normalization, added FMV, 07/24/03
597  Double_t xpmin = x.min(rangeName)/rtau ;
598  Double_t xpmax = x.max(rangeName)/rtau ;
599  Double_t c = sig/(root2*rtau) ;
600  Double_t umin = xpmin/(2*c) ;
601  Double_t umax = xpmax/(2*c) ;
602  Double_t result ;
603  if (_asympInt) {
604  result = 1.0 ;
605  } else {
606  result = 0.5*evalCerfInt(-fsign,rtau,fsign*umin,fsign*umax,c)/rtau ; //WVEFIX add 1/rtau
607  }
608 
609  if (_basisCode!=0 && basisSign==Both) result *= 2 ;
610  //cout << "Integral 1st form " << " result= " << result*ssfInt << endl;
611  return result*ssfInt ;
612  }
613 
614  Double_t omega = (basisType!=expBasis) ?((RooAbsReal*)basis().getParameter(2))->getVal() : 0 ;
615 
616  // *** 2nd form: unity, used for sinBasis and cosBasis with tau=0 (PDF is zero) ***
617  //if (tau==0&&omega!=0) {
618  if (tau==0) { // modified, FMV 07/24/03
619  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 2nd form" << endl ;
620  return 0. ;
621  }
622 
623  // *** 3rd form: Convolution with exp(-t/tau), used for expBasis and cosBasis(omega=0) ***
624  if (basisType==expBasis || (basisType==cosBasis && omega==0.)) {
625  //Double_t result = 2*tau ;
626  //if (basisSign==Both) result *= 2 ;
627  // finite+asymtotic normalization, added FMV, 07/24/03
628  Double_t result(0.);
629  if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,sig,rtau,fsign,rangeName);
630  if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,sig,rtau,fsign,rangeName);
631  //cout << "Integral 3rd form " << " result= " << result*ssfInt << endl;
632  return result*ssfInt ;
633  }
634 
635  // *** 4th form: Convolution with exp(-t/tau)*sin(omega*t), used for sinBasis(omega<>0,tau<>0) ***
636  Double_t wt = omega * tau ;
637  if (basisType==sinBasis) {
638  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName() << ") 4th form omega = "
639  << omega << ", tau = " << tau << endl ;
640  //cout << "sin integral" << endl;
641  Double_t result(0) ;
642  if (wt==0) return result ;
643  //if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega).imag() ;
644  //if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega).imag() ;
645  // finite+asymtotic normalization, added FMV, 07/24/03
646  if (basisSign!=Minus) result += -1*calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).imag();
647  if (basisSign!=Plus) result += -1*calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).imag();
648  //cout << "Integral 4th form " << " result= " << result*ssfInt << endl;
649  return result*ssfInt ;
650  }
651 
652  // *** 5th form: Convolution with exp(-t/tau)*cos(omega*t), used for cosBasis(omega<>0) ***
653  if (basisType==cosBasis) {
654  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
655  << ") 5th form omega = " << omega << ", tau = " << tau << endl ;
656  //cout << "cos integral" << endl;
657  Double_t result(0) ;
658  //if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega).real() ;
659  //if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega).real() ;
660  // finite+asymtotic normalization, added FMV, 07/24/03
661  if (basisSign!=Minus) result += calcSinConvNorm(+1,tau,omega,sig,rtau,fsign,rangeName).real();
662  if (basisSign!=Plus) result += calcSinConvNorm(-1,tau,omega,sig,rtau,fsign,rangeName).real();
663  //cout << "Integral 5th form " << " result= " << result*ssfInt << endl;
664  return result*ssfInt ;
665  }
666 
667  Double_t dgamma = ((basisType==coshBasis)||(basisType==sinhBasis))?((RooAbsReal*)basis().getParameter(2))->getVal():0 ;
668 
669  // *** 6th form: Convolution with exp(-t/tau)*sinh(dgamma*t/2), used for sinhBasis ***
670  if (basisType==sinhBasis) {
671  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
672  << ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
673  Double_t tau1 = 1/(1/tau-dgamma/2);
674  Double_t tau2 = 1/(1/tau+dgamma/2);
675  //cout << "sinh integral" << endl;
676  Double_t result(0) ;
677  //if (basisSign!=Minus) result += tau1-tau2 ;
678  //if (basisSign!=Plus) result += tau2-tau1 ;
679  // finite+asymtotic normalization, added FMV, 07/24/03
680  if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)-
681  calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
682  if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName)-
683  calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName));
684  //cout << "Integral 6th form " << " result= " << result*ssfInt << endl;
685  return result;
686  }
687 
688  // ** 7th form: Convolution with exp(-t/tau)*cosh(dgamma*t/2), used for coshBasis ***
689  if (basisType==coshBasis) {
690  if (verboseEval()>0) cout << "RooGExpModel::analyticalIntegral(" << GetName()
691  << ") 6th form dgamma = " << dgamma << ", tau = " << tau << endl ;
692  //cout << "cosh integral" << endl;
693  Double_t tau1 = 1/(1/tau-dgamma/2);
694  Double_t tau2 = 1/(1/tau+dgamma/2);
695  //Double_t result = (tau1+tau2) ;
696  //if (basisSign==Both) result *= 2 ;
697  // finite+asymtotic normalization, added FMV, 07/24/03
698  Double_t result(0);
699  if (basisSign!=Minus) result += 0.5*(calcSinConvNorm(+1,tau1,sig,rtau,fsign,rangeName)+
700  calcSinConvNorm(+1,tau2,sig,rtau,fsign,rangeName));
701  if (basisSign!=Plus) result += 0.5*(calcSinConvNorm(-1,tau1,sig,rtau,fsign,rangeName)+
702  calcSinConvNorm(-1,tau2,sig,rtau,fsign,rangeName));
703  //cout << "Integral 7th form " << " result= " << result*ssfInt << endl;
704  return result;
705 
706  }
707 
708  R__ASSERT(0) ;
709  return 1 ;
710 }
711 
712 
713 // modified FMV, 07/24/03. Finite+asymtotic normalization
714 
715 ////////////////////////////////////////////////////////////////////////////////
716 /// old code (asymptotic normalization only)
717 /// std::complex<Double_t> z(1/tau,sign*omega);
718 /// return z*2/(omega*omega+1/(tau*tau));
719 
720 std::complex<Double_t> RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t omega,
721  Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
722 {
723  static Double_t root2(sqrt(2.)) ;
724 
725  Double_t smin1= x.min(rangeName)/tau;
726  Double_t smax1= x.max(rangeName)/tau;
727  Double_t c1= sig/(root2*tau);
728  Double_t umin1= smin1/(2*c1);
729  Double_t umax1= smax1/(2*c1);
730  Double_t smin2= x.min(rangeName)/rtau;
731  Double_t smax2= x.max(rangeName)/rtau;
732  Double_t c2= sig/(root2*rtau);
733  Double_t umin2= smin2/(2*c2) ;
734  Double_t umax2= smax2/(2*c2) ;
735 
736  std::complex<Double_t> eins(1,0);
737  std::complex<Double_t> k(1/tau,sign*omega);
738  std::complex<Double_t> term1 = evalCerfInt(sign,-sign*omega*tau, tau, -sign*umin1, -sign*umax1, c1);
739  //std::complex<Double_t> term2 = evalCerfInt(-fsign,0., rtau, fsign*umin2, fsign*umax2, c2)*std::complex<Double_t>(fsign*sign,0);
740  std::complex<Double_t> term2 = std::complex<Double_t>(evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign,0);
741  return (term1+term2)/(eins + k*fsign*sign*rtau) ;
742 }
743 
744 
745 // added FMV, 08/17/03
746 
747 ////////////////////////////////////////////////////////////////////////////////
748 
749 Double_t RooGExpModel::calcSinConvNorm(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign, const char* rangeName) const
750 {
751  static Double_t root2(sqrt(2.)) ;
752 
753  Double_t smin1= x.min(rangeName)/tau;
754  Double_t smax1= x.max(rangeName)/tau;
755  Double_t c1= sig/(root2*tau);
756  Double_t umin1= smin1/(2*c1);
757  Double_t umax1= smax1/(2*c1);
758  Double_t smin2= x.min(rangeName)/rtau;
759  Double_t smax2= x.max(rangeName)/rtau;
760  Double_t c2= sig/(root2*rtau);
761  Double_t umin2= smin2/(2*c2) ;
762  Double_t umax2= smax2/(2*c2) ;
763 
764  Double_t eins(1);
765  Double_t k(1/tau);
766  Double_t term1 = evalCerfInt(sign, tau, -sign*umin1, -sign*umax1, c1);
767  Double_t term2 = evalCerfInt(-fsign, rtau, fsign*umin2, fsign*umax2, c2)*fsign*sign;
768 
769  // WVE Handle 0/0 numeric divergence
770  if (fabs(tau-rtau)<1e-10 && fabs(term1+term2)<1e-10) {
771  cout << "epsilon method" << endl ;
772  static Double_t epsilon = 1e-4 ;
773  return calcSinConvNorm(sign,tau+epsilon,sig,rtau-epsilon,fsign,rangeName) ;
774  }
775  return (term1+term2)/(eins + k*fsign*sign*rtau) ;
776 }
777 
778 
779 
780 // added FMV, 07/24/03
781 ////////////////////////////////////////////////////////////////////////////////
782 
783 std::complex<Double_t> RooGExpModel::evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
784 {
785  std::complex<Double_t> diff;
786  if (_asympInt) {
787  diff = std::complex<Double_t>(2,0) ;
788  } else {
789  diff = std::complex<Double_t>(sign,0.)*(evalCerf(wt,umin,c) - evalCerf(wt,umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
790  }
791  return std::complex<Double_t>(tau/(1.+wt*wt),0)*std::complex<Double_t>(1,wt)*diff;
792 }
793 // added FMV, 08/17/03. Modified FMV, 08/30/03
794 
795 ////////////////////////////////////////////////////////////////////////////////
796 
798 {
799  Double_t diff;
800  if (_asympInt) {
801  diff = 2. ;
802  } else {
803  if ((umin<-8 && umax>8)||(umax<-8 && umin>8)) {
804  // If integral is over >8 sigma, approximate with full integral
805  diff = 2. ;
806  } else {
807  diff = sign*(evalCerfRe(umin,c) - evalCerfRe(umax,c) + RooMath::erf(umin) - RooMath::erf(umax));
808  }
809  }
810  return tau*diff;
811 }
812 
813 
814 
815 ////////////////////////////////////////////////////////////////////////////////
816 /// use the approximation: erf(z) = exp(-z*z)/(sqrt(pi)*z)
817 /// to explicitly cancel the divergent exp(y*y) behaviour of
818 /// CWERF for z = x + i y with large negative y
819 
820 std::complex<Double_t> RooGExpModel::evalCerfApprox(Double_t swt, Double_t u, Double_t c)
821 {
822  static Double_t rootpi= sqrt(atan2(0.,-1.));
823  std::complex<Double_t> z(swt*c,u+c);
824  std::complex<Double_t> zc(u+c,-swt*c);
825  std::complex<Double_t> zsq= z*z;
826  std::complex<Double_t> v= -zsq - u*u;
827 
828  return std::exp(v)*(-std::exp(zsq)/(zc*rootpi) + 1.)*2.;
829 }
830 
831 
832 
833 ////////////////////////////////////////////////////////////////////////////////
834 
835 Int_t RooGExpModel::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
836 {
837  if (matchArgs(directVars,generateVars,x)) return 1 ;
838  return 0 ;
839 }
840 
841 
842 
843 ////////////////////////////////////////////////////////////////////////////////
844 
846 {
847  R__ASSERT(code==1) ;
848  Double_t xgen ;
849  while(1) {
851  Double_t xexp = RooRandom::uniform();
852  if (!_flip) xgen= xgau + (rlife*rsf)*log(xexp); // modified, FMV 08/13/03
853  else xgen= xgau - (rlife*rsf)*log(xexp);
854  if (xgen<x.max() && xgen>x.min()) {
855  x = xgen ;
856  return ;
857  }
858  }
859 }
860 
861 
862 
863 
864 
865 
866 
867 
868 
869 
870 
871 
872 
873 
874 
875 
876 
877 
878 
879 
880 
881 
882 
883 
884 
885 
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
XYZVector ans(TestRotation const &t, XYZVector const &v_in)
Double_t logErfC(Double_t x) const
Approximation of the log of the complex error function.
return c
virtual Double_t Gaus(Double_t mean=0, Double_t sigma=1)
Samples a random number from the standard Normal (Gaussian) Distribution with the given mean and sigm...
Definition: TRandom.cxx:235
return c1
Definition: legend1.C:41
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
#define R__ASSERT(e)
Definition: TError.h:98
RooRealVar & convVar() const
Return the convolution variable of the resolution model.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
TArc * a
Definition: textangle.C:12
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
virtual Int_t basisCode(const char *name) const
RooRealConstant provides static functions to create and keep track of RooRealVar constants.
static int verboseEval()
Return global level of verbosity for p.d.f. evaluations.
Definition: RooAbsPdf.cxx:2950
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:501
const RooFormulaVar & basis() const
double sqrt(double)
static std::complex< Double_t > evalCerfApprox(Double_t swt, Double_t u, Double_t c)
use the approximation: erf(z) = exp(-z*z)/(sqrt(pi)*z) to explicitly cancel the divergent exp(y*y) be...
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
friend class RooArgSet
Definition: RooAbsArg.h:469
const Double_t sigma
virtual Double_t evaluate() const
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:37
RooRealProxy rsf
Definition: RooGExpModel.h:119
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
SVector< double, 2 > v
Definition: Dict.h:5
RooRealProxy ssf
Definition: RooGExpModel.h:118
Double_t calcDecayConv(Double_t sign, Double_t tau, Double_t sig, Double_t rtau, Double_t fsign) const
static std::complex< double > erfc(const std::complex< double > z)
complex erfc function
Definition: RooMath.cxx:560
Bool_t _flatSFInt
Definition: RooGExpModel.h:122
REAL epsilon
Definition: triangle.c:617
static std::complex< Double_t > evalCerf(Double_t swt, Double_t u, Double_t c)
Definition: RooGExpModel.h:94
return c2
Definition: legend2.C:14
#define ClassImp(name)
Definition: Rtypes.h:279
static Double_t uniform(TRandom *generator=randomGenerator())
Return a number uniformly distributed from (0,1)
Definition: RooRandom.cxx:84
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Double_t evalCerfRe(Double_t u, Double_t c) const
Definition: RooGExpModel.h:102
Bool_t _asympInt
Definition: RooGExpModel.h:123
double atan2(double, double)
virtual ~RooGExpModel()
Destructor.
int type
Definition: TGX11.cxx:120
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
RooRealProxy rlife
Definition: RooGExpModel.h:117
std::complex< Double_t > evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
std::complex< Double_t > calcSinConvNorm(Double_t sign, Double_t tau, Double_t omega, Double_t sig, Double_t rtau, Double_t fsign, const char *rangeName) const
old code (asymptotic normalization only) std::complex<Double_t> z(1/tau,sign*omega); return z*2/(omeg...
Class RooGExpModel is a RooResolutionModel implementation that models a resolution function that is t...
Definition: RooGExpModel.h:27
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:584
double result[121]
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
double exp(double)
RooAbsArg * getParameter(const char *name) const
Definition: RooFormulaVar.h:39
return
Definition: HLFactory.cxx:514
char name[80]
Definition: TGX11.cxx:109
double log(double)
RooRealProxy sigma
Definition: RooGExpModel.h:116
std::complex< Double_t > calcSinConv(Double_t sign, Double_t sig, Double_t tau, Double_t omega, Double_t rtau, Double_t fsign) const