Logo ROOT  
Reference Guide
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 "BatchHelpers.h"
32#include "RooVDTHeaders.h"
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///cout << "exp(x=" << x << ",c=" << c << ")=" << exp(c*x) << endl ;
59
61 return exp(c*x);
62}
63
64////////////////////////////////////////////////////////////////////////////////
65
66Int_t RooExponential::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
67{
68 if (matchArgs(allVars,analVars,x)) return 1;
69 if (matchArgs(allVars,analVars,c)) return 2;
70 return 0 ;
71}
72
73////////////////////////////////////////////////////////////////////////////////
74
75Double_t RooExponential::analyticalIntegral(Int_t code, const char* rangeName) const
76{
77 assert(code == 1 || code ==2);
78
79 auto& constant = code == 1 ? c : x;
80 auto& integrand = code == 1 ? x : c;
81
82 if (constant == 0.0) {
83 return integrand.max(rangeName) - integrand.min(rangeName);
84 }
85
86 return (exp(constant*integrand.max(rangeName)) - exp(constant*integrand.min(rangeName)))
87 / constant;
88}
89
90
91namespace {
92
93template<class Tx, class Tc>
94void compute(size_t n, double* __restrict output, Tx x, Tc c) {
95
96 for (size_t i = 0; i < n; ++i) { //CHECK_VECTORISE
97 output[i] = _rf_fast_exp(x[i]*c[i]);
98 }
99}
100
101}
102
103////////////////////////////////////////////////////////////////////////////////
104/// Evaluate the exponential without normalising it on the given batch.
105/// \param[in] batchIndex Index of the batch to be computed.
106/// \param[in] batchSize Size of each batch. The last batch may be smaller.
107/// \return A span with the computed values.
108
109RooSpan<double> RooExponential::evaluateBatch(std::size_t begin, std::size_t batchSize) const {
110 using namespace BatchHelpers;
111 auto xData = x.getValBatch(begin, batchSize);
112 auto cData = c.getValBatch(begin, batchSize);
113 const bool batchX = !xData.empty();
114 const bool batchC = !cData.empty();
115
116 if (!batchX && !batchC) {
117 return {};
118 }
119 batchSize = findSize({ xData, cData });
120 auto output = _batchData.makeWritableBatchUnInit(begin, batchSize);
121
122 if (batchX && !batchC ) {
123 compute(batchSize, output.data(), xData, BracketAdapter<double>(c));
124 }
125 else if (!batchX && batchC ) {
126 compute(batchSize, output.data(), BracketAdapter<double>(x), cData);
127 }
128 else if (batchX && batchC ) {
129 compute(batchSize, output.data(), xData, cData);
130 }
131 return output;
132}
#define c(i)
Definition: RSha256.hxx:101
double _rf_fast_exp(double x)
VDT headers for RooFit.
Definition: RooVDTHeaders.h:47
#define ClassImp(name)
Definition: Rtypes.h:361
char name[80]
Definition: TGX11.cxx:109
double exp(double)
RooSpan< double > makeWritableBatchUnInit(std::size_t begin, std::size_t batchSize, const RooArgSet *const normSet=nullptr, Tag_t ownerTag=kUnspecified)
Make a batch and return a span pointing to the pdf-local memory.
Definition: BatchData.cxx:118
Little adapter that gives a bracket operator to types that don't have one.
Definition: BatchHelpers.h:58
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
BatchHelpers::BatchData _batchData
Definition: RooAbsReal.h:450
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Exponential PDF.
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.
RooSpan< double > evaluateBatch(std::size_t begin, std::size_t batchSize) const override
Evaluate the exponential without normalising it on the given batch.
RooRealProxy x
Double_t evaluate() const override
cout << "exp(x=" << x << ",c=" << c << ")=" << exp(c*x) << endl ;
A simple container to hold a batch of data values.
Definition: RooSpan.h:32
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
RooSpan< const typename T::value_type > getValBatch(std::size_t begin, std::size_t batchSize) const
Retrieve a batch of real or category data.
Double_t x[n]
Definition: legend1.C:17
const Int_t n
Definition: legend1.C:16
size_t findSize(std::vector< RooSpan< const double > > parameters)
This function returns the minimum size of the non-zero-sized batches.
static void output(int code)
Definition: gifencode.c:226