Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooUniform.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 RooUniform
18 \ingroup Roofit
19
20Flat p.d.f. in N dimensions
21**/
22
23#include "RooUniform.h"
24#include "RooBatchCompute.h"
25#include "RooAbsReal.h"
26#include "RooRealVar.h"
27#include "RooArgSet.h"
28
30
31////////////////////////////////////////////////////////////////////////////////
32
33RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
34 RooAbsPdf(name,title),
35 x("x","Observables",this,kTRUE,kFALSE)
36{
37 x.add(_x) ;
38}
39
40////////////////////////////////////////////////////////////////////////////////
41
42RooUniform::RooUniform(const RooUniform& other, const char* name) :
43 RooAbsPdf(other,name), x("x",this,other.x)
44{
45}
46
47////////////////////////////////////////////////////////////////////////////////
48
50{
51 return 1 ;
52}
53
54////////////////////////////////////////////////////////////////////////////////
55///Compute multiple values of the uniform distribution (effectively return a span with ones)
57{
58 size_t nEvents = 1;
59 for (auto elm : x) {
60 size_t nEventsCurrent = static_cast<const RooAbsReal*>(elm)->getValues(evalData).size();
61 if(nEventsCurrent != 1 && nEvents != 1 && nEventsCurrent != nEvents) {
62 auto errorMsg = std::string("RooUniform::evaluateSpan(): number of entries for input variables does not match")
63 + "in RooUniform with name \"" + GetName() + "\".";
64 coutE(FastEvaluations) << errorMsg << std::endl ;
65 throw std::runtime_error(errorMsg);
66 }
67 nEvents = std::max(nEvents, nEventsCurrent);
68 }
69 RooSpan<double> values = evalData.makeBatch(this, nEvents);
70 for (size_t i=0; i<nEvents; i++) {
71 values[i] = 1.0;
72 }
73 return values;
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Advertise analytical integral
78
79Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
80{
81 Int_t nx = x.getSize() ;
82 if (nx>31) {
83 // Warn that analytical integration is only provided for the first 31 observables
84 coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
85 << " observables, analytical integration is only implemented for the first 31 observables" << std::endl ;
86 nx=31 ;
87 }
88
89 Int_t code(0) ;
90 for (int i=0 ; i<x.getSize() ; i++) {
91 if (allVars.find(x.at(i)->GetName())) {
92 code |= (1<<i) ;
93 analVars.add(*allVars.find(x.at(i)->GetName())) ;
94 }
95 }
96 return code ;
97}
98
99////////////////////////////////////////////////////////////////////////////////
100/// Implement analytical integral
101
102Double_t RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
103{
104 Double_t ret(1) ;
105 for (int i=0 ; i<32 ; i++) {
106 if (code&(1<<i)) {
108 ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
109 }
110 }
111 return ret ;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Advertise internal generator
116
117Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
118{
119 Int_t nx = x.getSize() ;
120 if (nx>31) {
121 // Warn that analytical integration is only provided for the first 31 observables
122 coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
123 << " observables, internal integrator is only implemented for the first 31 observables" << std::endl ;
124 nx=31 ;
125 }
126
127 Int_t code(0) ;
128 for (int i=0 ; i<x.getSize() ; i++) {
129 if (directVars.find(x.at(i)->GetName())) {
130 code |= (1<<i) ;
131 generateVars.add(*directVars.find(x.at(i)->GetName())) ;
132 }
133 }
134 return code ;
135 return 0 ;
136}
137
138////////////////////////////////////////////////////////////////////////////////
139/// Implement internal generator
140
142{
143 // Fast-track handling of one-observable case
144 if (code==1) {
145 ((RooAbsRealLValue*)x.at(0))->randomize() ;
146 return ;
147 }
148
149 for (int i=0 ; i<32 ; i++) {
150 if (code&(1<<i)) {
152 var->randomize() ;
153 }
154 }
155}
#define coutW(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:92
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:91
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
Int_t getSize() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute batch of values for given input data, and normalise by integrating over the observables in no...
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
virtual Double_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual void randomize(const char *rangeName=0)
Set a new value sampled from a uniform distribution over the fit range.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:70
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Add element to non-owning set.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Reimplementation of standard RooArgList::add()
A simple container to hold a batch of data values.
Definition RooSpan.h:34
Flat p.d.f.
Definition RooUniform.h:24
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *=nullptr) const
Compute multiple values of the uniform distribution (effectively return a span with ones)
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Advertise analytical integral.
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implement analytical integral.
void generateEvent(Int_t code)
Implement internal generator.
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Advertise internal generator.
RooListProxy x
Definition RooUniform.h:40
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31
RooSpan< double > makeBatch(const RooAbsReal *owner, std::size_t size)
Create a writable batch.