Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
30
31#include "TMath.h"
32#include "Math/ProbFunc.h"
33
34
35////////////////////////////////////////////////////////////////////////////////
36
37RooLandau::RooLandau(const char *name, const char *title, RooAbsReal::Ref _x, RooAbsReal::Ref _mean, RooAbsReal::Ref _sigma) :
38 RooAbsPdf(name,title),
39 x("x","Dependent",this,_x),
40 mean("mean","Mean",this,_mean),
41 sigma("sigma","Width",this,_sigma)
42{
43 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal&>(_sigma)}, 0.0);
44}
45
46////////////////////////////////////////////////////////////////////////////////
47
50 x("x",this,other.x),
51 mean("mean",this,other.mean),
52 sigma("sigma",this,other.sigma)
53{
54}
55
56////////////////////////////////////////////////////////////////////////////////
57
62
63////////////////////////////////////////////////////////////////////////////////
64/// Compute multiple values of Landau distribution.
66{
68 {ctx.at(x), ctx.at(mean), ctx.at(sigma)});
69}
70
71////////////////////////////////////////////////////////////////////////////////
72
74{
75 return matchArgs(directVars,generateVars,x) ? 1 : 0;
76}
77
78Int_t RooLandau::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
79{
80 return matchArgs(allVars, analVars, x) ? 1 : 0;
81}
82
84{
85 // Don't do anything with "code". It can only be "1" anyway (see
86 // implementation of getAnalyticalIntegral).
87
88 const double meanVal = mean;
89 const double sigmaVal = sigma;
90
93 return sigmaVal * (a - b);
94}
95
96////////////////////////////////////////////////////////////////////////////////
97
99{
100 assert(1 == code); (void)code;
101 double xgen ;
102 while(true) {
104 if (xgen<x.max() && xgen>x.min()) {
105 x = xgen ;
106 break;
107 }
108 }
109 return;
110}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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:24
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:78
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition RooLandau.cxx:98
RooRealProxy x
Definition RooLandau.h:46
RooRealProxy sigma
Definition RooLandau.h:48
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition RooLandau.cxx:58
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Landau distribution.
Definition RooLandau.cxx:65
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:73
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:83
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
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.
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={})
double landau(double x, double mu, double sigma)
Definition MathFuncs.h:357
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.