Skip to content

Commit 52b25ef

Browse files
committed
FEAT: RooFit - Add 1D bounded constructor to RooUniform
Signed-off-by: JasMehta08 <[email protected]>
1 parent 2e2dfc5 commit 52b25ef

File tree

2 files changed

+182
-79
lines changed

2 files changed

+182
-79
lines changed

roofit/roofit/inc/RooUniform.h

Lines changed: 18 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -18,32 +18,36 @@
1818

1919
#include "RooAbsPdf.h"
2020
#include "RooListProxy.h"
21+
#include "RooRealProxy.h"
2122

2223
class RooRealVar;
24+
class RooAbsReal;
2325

2426
class RooUniform : public RooAbsPdf {
2527
public:
26-
RooUniform() {} ;
27-
RooUniform(const char *name, const char *title, const RooArgSet& _x);
28-
RooUniform(const RooUniform& other, const char* name=nullptr) ;
29-
TObject* clone(const char* newname=nullptr) const override { return new RooUniform(*this,newname); }
28+
RooUniform() {};
29+
RooUniform(const char *name, const char *title, const RooArgSet &_x);
30+
/// Constructor for a 1D uniform PDF with fittable bounds.
31+
RooUniform(const char *name, const char *title, RooAbsReal &x, RooAbsReal &x_low, RooAbsReal &x_up);
32+
RooUniform(const RooUniform &other, const char *name = nullptr);
33+
TObject *clone(const char *newname = nullptr) const override { return new RooUniform(*this, newname); }
3034

31-
Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
32-
double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
35+
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName = nullptr) const override;
36+
double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override;
3337

34-
Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK=true) const override;
35-
void generateEvent(Int_t code) override;
38+
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK = true) const override;
39+
void generateEvent(Int_t code) override;
3640

3741
protected:
42+
RooListProxy x; ///< List of observables for N-dimensional uniform PDF
43+
RooRealProxy x_single; ///< Single observable for 1D bounded mode
44+
RooRealProxy x_low; ///< Lower bound for 1D bounded mode
45+
RooRealProxy x_up; ///< Upper bound for 1D bounded mode
3846

39-
RooListProxy x ;
40-
41-
double evaluate() const override ;
42-
47+
double evaluate() const override;
4348

4449
private:
45-
46-
ClassDefOverride(RooUniform,1) // Flat PDF in N dimensions
50+
ClassDefOverride(RooUniform, 1) // Flat PDF in N dimensions
4751
};
4852

4953
#endif

roofit/roofit/src/RooUniform.cxx

Lines changed: 164 additions & 65 deletions
Original file line numberDiff line numberDiff line change
@@ -17,114 +17,213 @@
1717
/** \class RooUniform
1818
\ingroup Roofit
1919
20-
Flat p.d.f. in N dimensions
20+
Flat p.d.f. in N dimensions or in 1D with explicit, fittable bounds.
21+
22+
This class can be used in two ways:
23+
1. **Multi-dimensional (Legacy):** By providing a RooArgSet of observables in the constructor, it creates
24+
a PDF that is uniform over the full range of all observables. This is the original behavior.
25+
26+
2. **1D with Fittable Bounds:** By providing a single observable (`x`) and two RooAbsReal
27+
parameters for the lower (`x_low`) and upper (`x_up`) bounds, it creates a 1D PDF that is uniform
28+
only between those bounds and zero everywhere else.
29+
30+
Example of the 1D bounded mode:
31+
```
32+
RooRealVar x("x", "x", 0, 10);
33+
RooRealVar low("low", "low", 2, 0, 10);
34+
RooRealVar high("high", "high", 8, 0, 10);
35+
RooUniform bounded_uniform("bounded", "bounded", x, low, high);
36+
```
2137
**/
2238

2339
#include "RooArgSet.h"
2440
#include "RooRealVar.h"
2541
#include "RooUniform.h"
42+
#include "RooRandom.h"
43+
#include <algorithm>
2644

27-
45+
ClassImp(RooUniform);
2846

2947
////////////////////////////////////////////////////////////////////////////////
30-
31-
RooUniform::RooUniform(const char *name, const char *title, const RooArgSet& _x) :
32-
RooAbsPdf(name,title),
33-
x("x","Observables",this,true,false)
48+
/// Legacy constructor for an N-dimensional uniform PDF.
49+
50+
RooUniform::RooUniform(const char *name, const char *title, const RooArgSet &vars)
51+
: RooAbsPdf(name, title),
52+
x("x", "Observables", this, true, false),
53+
x_single("x_single", "", this),
54+
x_low("x_low", "", this),
55+
x_up("x_up", "", this)
3456
{
35-
x.add(_x) ;
57+
x.add(vars);
3658
}
3759

3860
////////////////////////////////////////////////////////////////////////////////
61+
/// Constructor for a 1D uniform PDF with fittable bounds.
62+
63+
RooUniform::RooUniform(const char *name, const char *title, RooAbsReal &var, RooAbsReal &low, RooAbsReal &up)
64+
: RooAbsPdf(name, title),
65+
x("x", "Observables", this, true, false),
66+
x_single("x_single", "Observable", this, var),
67+
x_low("x_low", "Lower Bound", this, low),
68+
x_up("x_up", "Upper Bound", this, up)
69+
{
70+
// Add the single observable to the list proxy for compatibility with legacy code
71+
x.add(var);
72+
}
3973

40-
RooUniform::RooUniform(const RooUniform& other, const char* name) :
41-
RooAbsPdf(other,name), x("x",this,other.x)
74+
////////////////////////////////////////////////////////////////////////////////
75+
/// Copy constructor
76+
77+
RooUniform::RooUniform(const RooUniform &other, const char *name)
78+
: RooAbsPdf(other, name),
79+
x("x", this, other.x),
80+
x_single("x_single", this, other.x_single),
81+
x_low("x_low", this, other.x_low),
82+
x_up("x_up", this, other.x_up)
4283
{
4384
}
4485

4586
////////////////////////////////////////////////////////////////////////////////
4687

4788
double RooUniform::evaluate() const
4889
{
49-
return 1 ;
90+
// If 1D bounded mode (check if proxies are valid)
91+
if (x_single.absArg() && x_low.absArg() && x_up.absArg()) {
92+
double low = x_low;
93+
double up = x_up;
94+
95+
if (low >= up)
96+
return 0.0; // Unphysical bounds
97+
98+
double val = x_single;
99+
if (val >= low && val <= up) {
100+
// Return the correctly normalized value for a uniform PDF
101+
return 1.0;
102+
} else {
103+
return 0.0;
104+
}
105+
}
106+
// uniform over observable range
107+
return 1.0;
50108
}
51109

52110
////////////////////////////////////////////////////////////////////////////////
53111
/// Advertise analytical integral
54112

55-
Int_t RooUniform::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
113+
Int_t RooUniform::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
56114
{
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 ;
115+
// 1D explicit bounds mode
116+
if (x_single.absArg() && x_low.absArg() && x_up.absArg()) {
117+
if (matchArgs(allVars, analVars, x_single))
118+
return 1;
119+
return 0;
120+
}
121+
// multi-dimensional mode
122+
Int_t nx = x.size();
123+
if (nx > 31) {
124+
coutW(Integration) << "RooUniform::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << x.size()
125+
<< " observables, analytical integration is only implemented for the first 31 observables"
126+
<< std::endl;
127+
nx = 31;
128+
}
129+
Int_t code(0);
130+
for (std::size_t i = 0; i < x.size(); i++) {
131+
if (allVars.find(x.at(i)->GetName())) {
132+
code |= (1 << i);
133+
analVars.add(*allVars.find(x.at(i)->GetName()));
134+
}
135+
}
136+
return code;
73137
}
74138

75139
////////////////////////////////////////////////////////////////////////////////
76140
/// Implement analytical integral
77141

78-
double RooUniform::analyticalIntegral(Int_t code, const char* rangeName) const
142+
double RooUniform::analyticalIntegral(Int_t code, const char *rangeName) const
79143
{
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 ;
144+
// 1D explicit bounds mode
145+
if (code == 1 && x_single.absArg() && x_low.absArg() && x_up.absArg()) {
146+
const RooAbsRealLValue &var = static_cast<const RooAbsRealLValue &>(x_single.arg());
147+
double low = x_low;
148+
double up = x_up;
149+
150+
if (low >= up)
151+
return 0.0;
152+
153+
double xmin = std::max(var.getMin(rangeName), low);
154+
double xmax = std::min(var.getMax(rangeName), up);
155+
156+
if (xmax > xmin)
157+
return (xmax - xmin);
158+
return 0.0;
159+
}
160+
161+
// multi-dimensional mode
162+
double ret(1);
163+
for (int i = 0; i < 32; i++) {
164+
if (code & (1 << i)) {
165+
RooAbsRealLValue *var = static_cast<RooAbsRealLValue *>(x.at(i));
166+
ret *= (var->getMax(rangeName) - var->getMin(rangeName));
167+
}
168+
}
169+
return ret;
88170
}
89171

90172
////////////////////////////////////////////////////////////////////////////////
91173
/// Advertise internal generator
92174

93-
Int_t RooUniform::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
175+
Int_t RooUniform::getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
94176
{
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 ;
177+
if (x_single.absArg() && x_low.absArg() && x_up.absArg()) {
178+
if (matchArgs(directVars, generateVars, x_single))
179+
return 2;
180+
return 0;
181+
}
182+
Int_t nx = x.size();
183+
if (nx > 31) {
184+
// Warn that analytical integration is only provided for the first 31 observables
185+
coutW(Integration) << "RooUniform::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << x.size()
186+
<< " observables, internal integrator is only implemented for the first 31 observables"
187+
<< std::endl;
188+
nx = 31;
189+
}
190+
191+
Int_t code(0);
192+
for (std::size_t i = 0; i < x.size(); i++) {
193+
if (directVars.find(x.at(i)->GetName())) {
194+
code |= (1 << i);
195+
generateVars.add(*directVars.find(x.at(i)->GetName()));
196+
}
197+
}
198+
return code;
111199
}
112200

113201
////////////////////////////////////////////////////////////////////////////////
114202
/// Implement internal generator
115203

116204
void RooUniform::generateEvent(Int_t code)
117205
{
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-
}
206+
if (code == 2) { // 1D Bounded case
207+
double low = x_low;
208+
double up = x_up;
209+
if (low < up) {
210+
static_cast<RooAbsRealLValue *>(x_single.absArg())->setVal(low + (up - low) * RooRandom::uniform());
211+
}
212+
return;
213+
}
214+
215+
// multi-dimensional case
216+
217+
// Fast-track handling of one-observable case
218+
if (code == 1) {
219+
(static_cast<RooAbsRealLValue *>(x.at(0)))->randomize();
220+
return;
221+
}
222+
223+
for (int i = 0; i < 32; i++) {
224+
if (code & (1 << i)) {
225+
RooAbsRealLValue *var = static_cast<RooAbsRealLValue *>(x.at(i));
226+
var->randomize();
227+
}
228+
}
130229
}

0 commit comments

Comments
 (0)