Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooExponential.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 RooExponential
18 \ingroup Roofit
19
20Exponential PDF. It computes
21\f[
22 \mathrm{RooExponential}(x, c) = \mathcal{N} \cdot \exp(c\cdot x),
23\f]
24where \f$ \mathcal{N} \f$ is a normalisation constant that depends on the
25range and values of the arguments.
26**/
27
28#include "RooExponential.h"
29
30#include "RooRealVar.h"
31#include "RooBatchCompute.h"
32
33
34#include <cmath>
35
36using namespace std;
37
39
40////////////////////////////////////////////////////////////////////////////////
41
42RooExponential::RooExponential(const char *name, const char *title,
43 RooAbsReal& _x, RooAbsReal& _c) :
44 RooAbsPdf(name, title),
45 x("x","Dependent",this,_x),
46 c("c","Exponent",this,_c)
47{
48}
49
50////////////////////////////////////////////////////////////////////////////////
51
53 RooAbsPdf(other, name), x("x",this,other.x), c("c",this,other.c)
54{
55}
56
57////////////////////////////////////////////////////////////////////////////////
58
60 return exp(c*x);
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Compute multiple values of Exponential distribution.
65void RooExponential::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
66{
68 dispatch->compute(stream, RooBatchCompute::Exponential, output, nEvents, {dataMap.at(x),dataMap.at(c)});
69}
70
71
72Int_t RooExponential::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
73{
74 if (matchArgs(allVars,analVars,x)) return 1;
75 if (matchArgs(allVars,analVars,c)) return 2;
76 return 0 ;
77}
78
79////////////////////////////////////////////////////////////////////////////////
80
81Double_t RooExponential::analyticalIntegral(Int_t code, const char* rangeName) const
82{
83 assert(code == 1 || code ==2);
84
85 auto& constant = code == 1 ? c : x;
86 auto& integrand = code == 1 ? x : c;
87
88 if (constant == 0.0) {
89 return integrand.max(rangeName) - integrand.min(rangeName);
90 }
91
92 return (exp(constant*integrand.max(rangeName)) - exp(constant*integrand.min(rangeName)))
93 / constant;
94}
#define c(i)
Definition RSha256.hxx:101
#define ClassImp(name)
Definition Rtypes.h:364
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_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:35
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
Exponential PDF.
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Exponential distribution.
RooRealProxy c
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy x
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
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...
static void output(int code)
Definition gifencode.c:226