Logo ROOT   6.14/05
Reference Guide
RooLognormal.cxx
Go to the documentation of this file.
1  /*****************************************************************************
2  * Project: RooFit *
3  * @(#)root/roofit:$Id$ *
4  * *
5  * RooFit Lognormal PDF *
6  * *
7  * Author: Gregory Schott and Stefan Schmitz *
8  * *
9  *****************************************************************************/
10 
11 /** \class RooLognormal
12  \ingroup Roofit
13 
14 RooFit Lognormal PDF. The two parameters are:
15  - m0: the median of the distribution
16  - k=exp(sigma): sigma is called the shape parameter in the TMath parameterization
17 
18 \f[ Lognormal(x,m_0,k) = \frac{e^{(-ln^2(x/m_0))/(2ln^2(k))}}{\sqrt{2\pi \cdot ln(k)\cdot x}} \f]
19 
20 The parameterization here is physics driven and differs from the ROOT::Math::lognormal_pdf(x,m,s,x0) with:
21  - m = log(m0)
22  - s = log(k)
23  - x0 = 0
24 **/
25 
26 #include "RooFit.h"
27 
28 #include "Riostream.h"
29 #include "Riostream.h"
30 #include <math.h>
31 
32 #include "RooLognormal.h"
33 #include "RooAbsReal.h"
34 #include "RooRealVar.h"
35 #include "RooRandom.h"
36 #include "RooMath.h"
37 #include "TMath.h"
38 
39 #include <Math/SpecFuncMathCore.h>
40 #include <Math/PdfFuncMathCore.h>
41 #include <Math/ProbFuncMathCore.h>
42 
43 #include "TError.h"
44 
45 using namespace std;
46 
48 
49 ////////////////////////////////////////////////////////////////////////////////
50 
51 RooLognormal::RooLognormal(const char *name, const char *title,
52  RooAbsReal& _x, RooAbsReal& _m0,
53  RooAbsReal& _k) :
54  RooAbsPdf(name,title),
55  x("x","Observable",this,_x),
56  m0("m0","m0",this,_m0),
57  k("k","k",this,_k)
58 {
59 }
60 
61 ////////////////////////////////////////////////////////////////////////////////
62 
63 RooLognormal::RooLognormal(const RooLognormal& other, const char* name) :
64  RooAbsPdf(other,name), x("x",this,other.x), m0("m0",this,other.m0),
65  k("k",this,other.k)
66 {
67 }
68 
69 ////////////////////////////////////////////////////////////////////////////////
70 /// ln(k)<1 would correspond to sigma < 0 in the parameterization
71 /// resulting by transforming a normal random variable in its
72 /// standard parameterization to a lognormal random variable
73 /// => treat ln(k) as -ln(k) for k<1
74 
76 {
77  Double_t xv = x;
79  Double_t ln_m0 = TMath::Log(m0);
80  Double_t x0 = 0;
81 
82  Double_t ret = ROOT::Math::lognormal_pdf(xv,ln_m0,ln_k,x0);
83  return ret ;
84 }
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 
88 Int_t RooLognormal::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
89 {
90  if (matchArgs(allVars,analVars,x)) return 1 ;
91  return 0 ;
92 }
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 
96 Double_t RooLognormal::analyticalIntegral(Int_t code, const char* rangeName) const
97 {
98  R__ASSERT(code==1) ;
99 
100  static const Double_t root2 = sqrt(2.) ;
101 
102  Double_t ln_k = TMath::Abs(TMath::Log(k));
103  Double_t ret = 0.5*( RooMath::erf( TMath::Log(x.max(rangeName)/m0)/(root2*ln_k) ) - RooMath::erf( TMath::Log(x.min(rangeName)/m0)/(root2*ln_k) ) ) ;
104 
105  return ret ;
106 }
107 
108 ////////////////////////////////////////////////////////////////////////////////
109 
110 Int_t RooLognormal::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
111 {
112  if (matchArgs(directVars,generateVars,x)) return 1 ;
113  return 0 ;
114 }
115 
116 ////////////////////////////////////////////////////////////////////////////////
117 
119 {
120  R__ASSERT(code==1) ;
121 
122  Double_t xgen ;
123  while(1) {
125  if (xgen<=x.max() && xgen>=x.min()) {
126  x = xgen ;
127  break;
128  }
129  }
130 
131  return;
132 }
double lognormal_pdf(double x, double m, double s, double x0=0)
Probability density function of the lognormal distribution.
Double_t Log(Double_t x)
Definition: TMath.h:759
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
#define R__ASSERT(e)
Definition: TError.h:96
Double_t Gaus(Double_t x, Double_t mean, Double_t sigma)
Gauss.
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
RooRealProxy k
Definition: RooLognormal.h:38
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:108
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...
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
RooRealProxy x
Definition: RooLognormal.h:36
RooFit Lognormal PDF.
Definition: RooLognormal.h:19
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Double_t evaluate() const
ln(k)<1 would correspond to sigma < 0 in the parameterization resulting by transforming a normal rand...
RooRealProxy m0
Definition: RooLognormal.h:37
Double_t Exp(Double_t x)
Definition: TMath.h:726
#define ClassImp(name)
Definition: Rtypes.h:359
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 max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:580
char name[80]
Definition: TGX11.cxx:109