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 "RooFit.h"
27#include "RooRandom.h"
28#include "RooBatchCompute.h"
29
30#include "TMath.h"
31
33
34////////////////////////////////////////////////////////////////////////////////
35
36RooLandau::RooLandau(const char *name, const char *title, RooAbsReal& _x, RooAbsReal& _mean, RooAbsReal& _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, {&_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
58{
59 return TMath::Landau(x, mean, sigma);
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Compute multiple values of Landau distribution.
65 return RooBatchCompute::dispatch->computeLandau(this, evalData, x->getValues(evalData, normSet), mean->getValues(evalData, normSet), sigma->getValues(evalData, normSet));
66}
67
68////////////////////////////////////////////////////////////////////////////////
69
70Int_t RooLandau::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
71{
72 if (matchArgs(directVars,generateVars,x)) return 1 ;
73 return 0 ;
74}
75
76////////////////////////////////////////////////////////////////////////////////
77
79{
80 assert(1 == code); (void)code;
81 Double_t xgen ;
82 while(1) {
84 if (xgen<x.max() && xgen>x.min()) {
85 x = xgen ;
86 break;
87 }
88 }
89 return;
90}
double Double_t
Definition RtypesCore.h:59
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
typedef void((*Func_t)())
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
Bool_t 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:29
virtual RooSpan< double > computeLandau(const RooAbsReal *, RunContext &, RooSpan< const double > x, RooSpan< const double > mean, RooSpan< const double > sigma)=0
Landau distribution p.d.f.
Definition RooLandau.h:24
void generateEvent(Int_t code)
Interface for generation of an event using the algorithm corresponding to the specified code.
Definition RooLandau.cxx:78
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition RooLandau.cxx:57
RooRealProxy x
Definition RooLandau.h:37
RooRealProxy sigma
Definition RooLandau.h:39
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...
Definition RooLandau.cxx:70
RooRealProxy mean
Definition RooLandau.h:38
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of Landau distribution.
Definition RooLandau.cxx:64
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:53
A simple container to hold a batch of data values.
Definition RooSpan.h:34
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
R__EXTERN RooBatchComputeInterface * dispatch
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 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:469
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31