Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooLandau.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/** \class RooLandau
18 \ingroup Roofit
19
20Landau distribution p.d.f
21\image html RF_Landau.png "PDF of the Landau distribution."
22**/
23
24#include "RooLandau.h"
25#include "RooHelpers.h"
26#include "RooRandom.h"
27#include "RooBatchCompute.h"
28
29#include "TMath.h"
30#include "Math/ProbFunc.h"
31
33
34////////////////////////////////////////////////////////////////////////////////
35
36RooLandau::RooLandau(const char *name, const char *title, RooAbsReal::Ref _x, RooAbsReal::Ref _mean, RooAbsReal::Ref _sigma) :
37 RooAbsPdf(name,title),
38 x("x","Dependent",this,_x),
39 mean("mean","Mean",this,_mean),
40 sigma("sigma","Width",this,_sigma)
41{
42 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal&>(_sigma)}, 0.0);
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
47RooLandau::RooLandau(const RooLandau& other, const char* name) :
48 RooAbsPdf(other,name),
49 x("x",this,other.x),
50 mean("mean",this,other.mean),
51 sigma("sigma",this,other.sigma)
52{
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
57double RooLandau::evaluate() const
58{
59 return TMath::Landau(x, mean, sigma);
60}
61
62////////////////////////////////////////////////////////////////////////////////
63
65{
66 ctx.addResult(this, ctx.buildCall("TMath::Landau", x, mean, sigma));
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Compute multiple values of Landau distribution.
72{
74 {ctx.at(x), ctx.at(mean), ctx.at(sigma)});
75}
76
77////////////////////////////////////////////////////////////////////////////////
78
79Int_t RooLandau::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
80{
81 if (matchArgs(directVars,generateVars,x)) return 1 ;
82 return 0;
83}
84
85Int_t RooLandau::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
86{
87 if (matchArgs(allVars, analVars, x))
88 return 1;
89 return 0;
90}
91
92Double_t RooLandau::analyticalIntegral(Int_t /*code*/, const char *rangeName) const
93{
94 // Don't do anything with "code". It can only be "1" anyway (see
95 // implementation of getAnalyticalIntegral).
96
97 const double meanVal = mean;
98 const double sigmaVal = sigma;
99
100 const double a = ROOT::Math::landau_cdf(x.max(rangeName), sigmaVal, meanVal);
101 const double b = ROOT::Math::landau_cdf(x.min(rangeName), sigmaVal, meanVal);
102 return sigmaVal * (a - b);
103}
104
105////////////////////////////////////////////////////////////////////////////////
106
107std::string RooLandau::buildCallToAnalyticIntegral(Int_t /*code*/, const char *rangeName,
109{
110 // Don't do anything with "code". It can only be "1" anyway (see
111 // implementation of getAnalyticalIntegral).
112 const std::string a = ctx.buildCall("ROOT::Math::landau_cdf", x.max(rangeName), sigma, mean);
113 const std::string b = ctx.buildCall("ROOT::Math::landau_cdf", x.min(rangeName), sigma, mean);
114 return ctx.getResult(sigma) + " * " + "(" + a + " - " + b + ")";
115}
116
117////////////////////////////////////////////////////////////////////////////////
118
120{
121 assert(1 == code); (void)code;
122 double xgen ;
123 while(true) {
125 if (xgen<x.max() && xgen>x.min()) {
126 x = xgen ;
127 break;
128 }
129 }
130 return;
131}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:68
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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
A class to maintain the context for squashing of RooFit models into code.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
std::string const & getResult(RooAbsArg const &arg)
Gets the result for the given node using the node name.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
Landau distribution p.d.f.
Definition RooLandau.h:24
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition RooLandau.cxx:85
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
RooRealProxy x
Definition RooLandau.h:46
RooRealProxy sigma
Definition RooLandau.h:48
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
Definition RooLandau.cxx:64
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition RooLandau.cxx:57
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Landau distribution.
Definition RooLandau.cxx:71
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...
Definition RooLandau.cxx:79
RooRealProxy mean
Definition RooLandau.h:47
double analyticalIntegral(Int_t code, const char *rangeName) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition RooLandau.cxx:92
std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:48
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.
virtual Double_t Landau(Double_t mean=0, Double_t sigma=1)
Generate a random number following a Landau distribution with location parameter mu and scale paramet...
Definition TRandom.cxx:381
double landau_cdf(double x, double xi=1, double x0=0)
Cumulative distribution function of the Landau distribution (lower tail).
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
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 Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition TMath.cxx:492