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 "RooArgSet.h"
24#include "RooRealVar.h"
25#include "RooUniform.h"
26
27
29
30////////////////////////////////////////////////////////////////////////////////
31
32RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
33 RooAbsPdf(name,title),
34 x("x","Observables",this,true,false)
35{
36 x.add(_x) ;
37}
38
39////////////////////////////////////////////////////////////////////////////////
40
41RooUniform::RooUniform(const RooUniform& other, const char* name) :
42 RooAbsPdf(other,name), x("x",this,other.x)
43{
44}
45
46////////////////////////////////////////////////////////////////////////////////
47
49{
50 return 1 ;
51}
52
53////////////////////////////////////////////////////////////////////////////////
54/// Advertise analytical integral
55
56Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
57{
58 Int_t nx = x.size() ;
59 if (nx>31) {
60 // Warn that analytical integration is only provided for the first 31 observables
61 coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.size()
62 << " observables, analytical integration is only implemented for the first 31 observables" << std::endl ;
63 nx=31 ;
64 }
65
66 Int_t code(0) ;
67 for (std::size_t i=0 ; i<x.size() ; i++) {
68 if (allVars.find(x.at(i)->GetName())) {
69 code |= (1<<i) ;
70 analVars.add(*allVars.find(x.at(i)->GetName())) ;
71 }
72 }
73 return code ;
74}
75
76////////////////////////////////////////////////////////////////////////////////
77/// Implement analytical integral
78
79double RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
80{
81 double ret(1) ;
82 for (int i=0 ; i<32 ; i++) {
83 if (code&(1<<i)) {
84 RooAbsRealLValue* var = static_cast<RooAbsRealLValue*>(x.at(i)) ;
85 ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
86 }
87 }
88 return ret ;
89}
90
91////////////////////////////////////////////////////////////////////////////////
92/// Advertise internal generator
93
94Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
95{
96 Int_t nx = x.size() ;
97 if (nx>31) {
98 // Warn that analytical integration is only provided for the first 31 observables
99 coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.size()
100 << " observables, internal integrator is only implemented for the first 31 observables" << std::endl ;
101 nx=31 ;
102 }
103
104 Int_t code(0) ;
105 for (std::size_t i=0 ; i<x.size() ; i++) {
106 if (directVars.find(x.at(i)->GetName())) {
107 code |= (1<<i) ;
108 generateVars.add(*directVars.find(x.at(i)->GetName())) ;
109 }
110 }
111 return code ;
112}
113
114////////////////////////////////////////////////////////////////////////////////
115/// Implement internal generator
116
118{
119 // Fast-track handling of one-observable case
120 if (code==1) {
121 (static_cast<RooAbsRealLValue*>(x.at(0)))->randomize() ;
122 return ;
123 }
124
125 for (int i=0 ; i<32 ; i++) {
126 if (code&(1<<i)) {
127 RooAbsRealLValue* var = static_cast<RooAbsRealLValue*>(x.at(i)) ;
128 var->randomize() ;
129 }
130 }
131}
#define coutW(a)
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
void randomize(const char *rangeName=nullptr) override
Set a new value sampled from a uniform distribution over the fit range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
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:55
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Flat p.d.f.
Definition RooUniform.h:24
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Advertise analytical integral.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implement analytical integral.
RooListProxy x
Definition RooUniform.h:39
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Advertise internal generator.
void generateEvent(Int_t code) override
Implement internal generator.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17