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 "RooAbsReal.h"
24#include "RooArgSet.h"
25#include "RooRealVar.h"
26#include "RooUniform.h"
27#include "RunContext.h"
28
29
31
32////////////////////////////////////////////////////////////////////////////////
33
34RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
35 RooAbsPdf(name,title),
36 x("x","Observables",this,kTRUE,kFALSE)
37{
38 x.add(_x) ;
39}
40
41////////////////////////////////////////////////////////////////////////////////
42
43RooUniform::RooUniform(const RooUniform& other, const char* name) :
44 RooAbsPdf(other,name), x("x",this,other.x)
45{
46}
47
48////////////////////////////////////////////////////////////////////////////////
49
51{
52 return 1 ;
53}
54
55////////////////////////////////////////////////////////////////////////////////
56///Compute multiple values of the uniform distribution (effectively return a span with ones)
58{
59 size_t nEvents = 1;
60 for (auto elm : x) {
61 size_t nEventsCurrent = static_cast<const RooAbsReal*>(elm)->getValues(evalData).size();
62 if(nEventsCurrent != 1 && nEvents != 1 && nEventsCurrent != nEvents) {
63 auto errorMsg = std::string("RooUniform::evaluateSpan(): number of entries for input variables does not match")
64 + "in RooUniform with name \"" + GetName() + "\".";
65 coutE(FastEvaluations) << errorMsg << std::endl ;
66 throw std::runtime_error(errorMsg);
67 }
68 nEvents = std::max(nEvents, nEventsCurrent);
69 }
70 RooSpan<double> values = evalData.makeBatch(this, nEvents);
71 for (size_t i=0; i<nEvents; i++) {
72 values[i] = 1.0;
73 }
74 return values;
75}
76
77////////////////////////////////////////////////////////////////////////////////
78/// Advertise analytical integral
79
80Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
81{
82 Int_t nx = x.getSize() ;
83 if (nx>31) {
84 // Warn that analytical integration is only provided for the first 31 observables
85 coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
86 << " observables, analytical integration is only implemented for the first 31 observables" << std::endl ;
87 nx=31 ;
88 }
89
90 Int_t code(0) ;
91 for (int i=0 ; i<x.getSize() ; i++) {
92 if (allVars.find(x.at(i)->GetName())) {
93 code |= (1<<i) ;
94 analVars.add(*allVars.find(x.at(i)->GetName())) ;
95 }
96 }
97 return code ;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101/// Implement analytical integral
102
103Double_t RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
104{
105 Double_t ret(1) ;
106 for (int i=0 ; i<32 ; i++) {
107 if (code&(1<<i)) {
109 ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
110 }
111 }
112 return ret ;
113}
114
115////////////////////////////////////////////////////////////////////////////////
116/// Advertise internal generator
117
118Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
119{
120 Int_t nx = x.getSize() ;
121 if (nx>31) {
122 // Warn that analytical integration is only provided for the first 31 observables
123 coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.getSize()
124 << " observables, internal integrator is only implemented for the first 31 observables" << std::endl ;
125 nx=31 ;
126 }
127
128 Int_t code(0) ;
129 for (int i=0 ; i<x.getSize() ; i++) {
130 if (directVars.find(x.at(i)->GetName())) {
131 code |= (1<<i) ;
132 generateVars.add(*directVars.find(x.at(i)->GetName())) ;
133 }
134 }
135 return code ;
136 return 0 ;
137}
138
139////////////////////////////////////////////////////////////////////////////////
140/// Implement internal generator
141
143{
144 // Fast-track handling of one-observable case
145 if (code==1) {
146 ((RooAbsRealLValue*)x.at(0))->randomize() ;
147 return ;
148 }
149
150 for (int i=0 ; i<32 ; i++) {
151 if (code&(1<<i)) {
153 var->randomize() ;
154 }
155 }
156}
#define coutW(a)
#define coutE(a)
const Bool_t kFALSE
Definition RtypesCore.h:101
double Double_t
Definition RtypesCore.h:59
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
Int_t getSize() const
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add the specified argument to list.
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:64
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
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:32
RooSpan< double > makeBatch(const RooAbsArg *owner, std::size_t size)
Create a writable batch.