Logo ROOT  
Reference Guide
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
32
33////////////////////////////////////////////////////////////////////////////////
34
35RooLandau::RooLandau(const char *name, const char *title, RooAbsReal& _x, RooAbsReal& _mean, RooAbsReal& _sigma) :
36 RooAbsPdf(name,title),
37 x("x","Dependent",this,_x),
38 mean("mean","Mean",this,_mean),
39 sigma("sigma","Width",this,_sigma)
40{
41 RooHelpers::checkRangeOfParameters(this, {&_sigma}, 0.0);
42}
43
44////////////////////////////////////////////////////////////////////////////////
45
46RooLandau::RooLandau(const RooLandau& other, const char* name) :
47 RooAbsPdf(other,name),
48 x("x",this,other.x),
49 mean("mean",this,other.mean),
50 sigma("sigma",this,other.sigma)
51{
52}
53
54////////////////////////////////////////////////////////////////////////////////
55
56double RooLandau::evaluate() const
57{
58 return TMath::Landau(x, mean, sigma);
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Compute multiple values of Landau distribution.
63void RooLandau::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
64{
66 dispatch->compute(stream, RooBatchCompute::Landau, output, nEvents,
67 {dataMap.at(x), dataMap.at(mean), dataMap.at(sigma)});
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
72Int_t RooLandau::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
73{
74 if (matchArgs(directVars,generateVars,x)) return 1 ;
75 return 0 ;
76}
77
78////////////////////////////////////////////////////////////////////////////////
79
81{
82 assert(1 == code); (void)code;
83 double xgen ;
84 while(1) {
86 if (xgen<x.max() && xgen>x.min()) {
87 x = xgen ;
88 break;
89 }
90 }
91 return;
92}
#define ClassImp(name)
Definition: Rtypes.h:375
char name[80]
Definition: TGX11.cxx:110
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:64
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:57
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition: DataMap.h:88
Landau distribution p.d.f.
Definition: RooLandau.h:24
RooLandau()
Definition: RooLandau.h:26
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition: RooLandau.cxx:80
RooRealProxy x
Definition: RooLandau.h:37
RooRealProxy sigma
Definition: RooLandau.h:39
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Landau distribution.
Definition: RooLandau.cxx:63
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooLandau.cxx:56
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:72
RooRealProxy mean
Definition: RooLandau.h:38
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:51
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper 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:380
const Double_t sigma
Double_t x[n]
Definition: legend1.C:17
void(off) SmallVectorTemplateBase< T
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.
Definition: RooHelpers.cxx:119
Double_t Landau(Double_t x, Double_t mpv=0, Double_t sigma=1, Bool_t norm=kFALSE)
The LANDAU function.
Definition: TMath.cxx:474
static void output()