Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
28
29////////////////////////////////////////////////////////////////////////////////
30
31RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
32 RooAbsPdf(name,title),
33 x("x","Observables",this,true,false)
34{
35 x.add(_x) ;
36}
37
38////////////////////////////////////////////////////////////////////////////////
39
42{
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
48{
49 return 1 ;
50}
51
52////////////////////////////////////////////////////////////////////////////////
53/// Advertise analytical integral
54
55Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
56{
57 Int_t nx = x.size() ;
58 if (nx>31) {
59 // Warn that analytical integration is only provided for the first 31 observables
60 coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.size()
61 << " observables, analytical integration is only implemented for the first 31 observables" << std::endl ;
62 nx=31 ;
63 }
64
65 Int_t code(0) ;
66 for (std::size_t i=0 ; i<x.size() ; i++) {
67 if (allVars.find(x.at(i)->GetName())) {
68 code |= (1<<i) ;
69 analVars.add(*allVars.find(x.at(i)->GetName())) ;
70 }
71 }
72 return code ;
73}
74
75////////////////////////////////////////////////////////////////////////////////
76/// Implement analytical integral
77
78double RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
79{
80 double ret(1) ;
81 for (int i=0 ; i<32 ; i++) {
82 if (code&(1<<i)) {
83 RooAbsRealLValue* var = static_cast<RooAbsRealLValue*>(x.at(i)) ;
84 ret *= (var->getMax(rangeName) - var->getMin(rangeName)) ;
85 }
86 }
87 return ret ;
88}
89
90////////////////////////////////////////////////////////////////////////////////
91/// Advertise internal generator
92
94{
95 Int_t nx = x.size() ;
96 if (nx>31) {
97 // Warn that analytical integration is only provided for the first 31 observables
98 coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.size()
99 << " observables, internal integrator is only implemented for the first 31 observables" << std::endl ;
100 nx=31 ;
101 }
102
103 Int_t code(0) ;
104 for (std::size_t i=0 ; i<x.size() ; i++) {
105 if (directVars.find(x.at(i)->GetName())) {
106 code |= (1<<i) ;
107 generateVars.add(*directVars.find(x.at(i)->GetName())) ;
108 }
109 }
110 return code ;
111}
112
113////////////////////////////////////////////////////////////////////////////////
114/// Implement internal generator
115
117{
118 // Fast-track handling of one-observable case
119 if (code==1) {
120 (static_cast<RooAbsRealLValue*>(x.at(0)))->randomize() ;
121 return ;
122 }
123
124 for (int i=0 ; i<32 ; i++) {
125 if (code&(1<<i)) {
126 RooAbsRealLValue* var = static_cast<RooAbsRealLValue*>(x.at(i)) ;
127 var->randomize() ;
128 }
129 }
130}
#define coutW(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
char name[80]
Definition TGX11.cxx:110
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:24
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:49
Double_t x[n]
Definition legend1.C:17