Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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
14RooFit 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[
19 \mathrm{RooLognormal}(x \, | \, m_0, k) = \frac{1}{\sqrt{2\pi \cdot \ln(k) \cdot x}} \cdot \exp\left( \frac{-\ln^2(\frac{x}{m_0})}{2 \ln^2(k)} \right)
20\f]
21
22The parameterization here is physics driven and differs from the ROOT::Math::lognormal_pdf() in `x,m,s,x0` with:
23 - `m = log(m0)`
24 - `s = log(k)`
25 - `x0 = 0`
26**/
27
28#include "RooLognormal.h"
29#include "RooRandom.h"
30#include "RooMath.h"
31#include "RooHelpers.h"
32#include "RooBatchCompute.h"
33
35
37
38////////////////////////////////////////////////////////////////////////////////
39
40RooLognormal::RooLognormal(const char *name, const char *title,
41 RooAbsReal& _x, RooAbsReal& _m0,
42 RooAbsReal& _k) :
43 RooAbsPdf(name,title),
44 x("x","Observable",this,_x),
45 m0("m0","m0",this,_m0),
46 k("k","k",this,_k)
47{
48 RooHelpers::checkRangeOfParameters(this, {&_x, &_m0, &_k}, 0.);
49
50 auto par = dynamic_cast<const RooAbsRealLValue*>(&_k);
51 if (par && par->getMin()<=1 && par->getMax()>=1 ) {
52 coutE(InputArguments) << "The parameter '" << par->GetName() << "' with range [" << par->getMin("") << ", "
53 << par->getMax() << "] of the " << this->ClassName() << " '" << this->GetName()
54 << "' can reach the unsafe value 1.0 " << ". Advise to limit its range." << std::endl;
55 }
56}
57
58////////////////////////////////////////////////////////////////////////////////
59
60RooLognormal::RooLognormal(const RooLognormal& other, const char* name) :
61 RooAbsPdf(other,name), x("x",this,other.x), m0("m0",this,other.m0),
62 k("k",this,other.k)
63{
64}
65
66////////////////////////////////////////////////////////////////////////////////
67/// ln(k)<1 would correspond to sigma < 0 in the parameterization
68/// resulting by transforming a normal random variable in its
69/// standard parameterization to a lognormal random variable
70/// => treat ln(k) as -ln(k) for k<1
71
73{
74 double ln_k = TMath::Abs(TMath::Log(k));
75 double ln_m0 = TMath::Log(m0);
76
77 double ret = ROOT::Math::lognormal_pdf(x,ln_m0,ln_k);
78 return ret ;
79}
80
81////////////////////////////////////////////////////////////////////////////////
82/// Compute multiple values of Lognormal distribution.
83void RooLognormal::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
84{
86 dispatch->compute(stream, RooBatchCompute::Lognormal, output, nEvents,
87 {dataMap.at(x), dataMap.at(m0), dataMap.at(k)});
88}
89
90////////////////////////////////////////////////////////////////////////////////
91
92Int_t RooLognormal::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
93{
94 if (matchArgs(allVars,analVars,x)) return 1 ;
95 return 0 ;
96}
97
98////////////////////////////////////////////////////////////////////////////////
99
100double RooLognormal::analyticalIntegral(Int_t code, const char* rangeName) const
101{
102 R__ASSERT(code==1) ;
103
104 static const double root2 = sqrt(2.) ;
105
106 double ln_k = TMath::Abs(TMath::Log(k));
107 double 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) ) ) ;
108
109 return ret ;
110}
111
112////////////////////////////////////////////////////////////////////////////////
113
114Int_t RooLognormal::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
115{
116 if (matchArgs(directVars,generateVars,x)) return 1 ;
117 return 0 ;
118}
119
120////////////////////////////////////////////////////////////////////////////////
121
123{
124 R__ASSERT(code==1) ;
125
126 double xgen ;
127 while(1) {
129 if (xgen<=x.max() && xgen>=x.min()) {
130 x = xgen ;
131 break;
132 }
133 }
134
135 return;
136}
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:117
char name[80]
Definition TGX11.cxx:110
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, ArgVector &)=0
RooSpan< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:21
RooFit Lognormal PDF.
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
RooRealProxy k
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Lognormal distribution.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy x
RooRealProxy m0
double evaluate() const override
ln(k)<1 would correspond to sigma < 0 in the parameterization resulting by transforming a normal rand...
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition RooMath.cxx:59
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:51
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:207
double lognormal_pdf(double x, double m, double s, double x0=0)
Probability density function of the lognormal distribution.
Double_t x[n]
Definition legend1.C:17
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string const &extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.
Double_t Exp(Double_t x)
Returns the base-e exponential function of x, which is e raised to the power x.
Definition TMath.h:707
Double_t Log(Double_t x)
Returns the natural logarithm of x.
Definition TMath.h:754
Short_t Abs(Short_t d)
Returns the absolute value of parameter Short_t d.
Definition TMathBase.h:123
static void output()