Skip to content

Commit 601ef3e

Browse files
committed
[RF] Update rf604_constraints to follow global observables convention
The RooFit convention for constraint terms is to have the constrained parameter also as the mean parameter of the constraint Gaussian, while a so-called global observable (a constant) is put in place of the pdf observable.
1 parent 997e03f commit 601ef3e

File tree

2 files changed

+12
-13
lines changed

2 files changed

+12
-13
lines changed

tutorials/roofit/rf604_constraints.C

Lines changed: 5 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -50,7 +50,7 @@ void rf604_constraints()
5050
// -----------------------------------------
5151

5252
// Construct Gaussian constraint pdf on parameter f at 0.8 with resolution of 0.1
53-
RooGaussian fconstraint("fconstraint", "fconstraint", f, RooConst(0.8), RooConst(0.2));
53+
RooGaussian fconstraint("fconstraint", "fconstraint", RooConst(0.8), f, RooConst(0.2));
5454

5555
// M E T H O D 1 - A d d i n t e r n a l c o n s t r a i n t t o m o d e l
5656
// -------------------------------------------------------------------------------------
@@ -62,19 +62,19 @@ void rf604_constraints()
6262
RooProdPdf modelc("modelc", "model with constraint", RooArgSet(model, fconstraint));
6363

6464
// Fit model (without use of constraint term)
65-
RooFitResult *r1 = model.fitTo(*d, Save());
65+
RooFitResult *r1 = model.fitTo(*d, Save(), PrintLevel(-1));
6666

6767
// Fit modelc with constraint term on parameter f
68-
RooFitResult *r2 = modelc.fitTo(*d, Constrain(f), Save());
68+
RooFitResult *r2 = modelc.fitTo(*d, Constrain(f), Save(), PrintLevel(-1));
6969

7070
// M E T H O D 2 - S p e c i f y e x t e r n a l c o n s t r a i n t w h e n f i t t i n g
7171
// -------------------------------------------------------------------------------------------------------
7272

7373
// Construct another Gaussian constraint pdf on parameter f at 0.2 with resolution of 0.1
74-
RooGaussian fconstext("fconstext", "fconstext", f, RooConst(0.2), RooConst(0.1));
74+
RooGaussian fconstext("fconstext", "fconstext", RooConst(0.2), f, RooConst(0.1));
7575

7676
// Fit with external constraint
77-
RooFitResult *r3 = model.fitTo(*d, ExternalConstraints(fconstext), Save());
77+
RooFitResult *r3 = model.fitTo(*d, ExternalConstraints(fconstext), Save(), PrintLevel(-1));
7878

7979
// Print the fit results
8080
cout << "fit result without constraint (data generated at f=0.5)" << endl;

tutorials/roofit/rf604_constraints.py

Lines changed: 7 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -37,33 +37,32 @@
3737

3838
# Construct Gaussian constraint pdf on parameter f at 0.8 with
3939
# resolution of 0.1
40-
fconstraint = ROOT.RooGaussian("fconstraint", "fconstraint", f, ROOT.RooFit.RooConst(0.8), ROOT.RooFit.RooConst(0.1))
40+
fconstraint = ROOT.RooGaussian("fconstraint", "fconstraint", ROOT.RooFit.RooConst(0.8), f, ROOT.RooFit.RooConst(0.1))
4141

4242
# Method 1 - add internal constraint to model
4343
# -------------------------------------------------------------------------------------
4444

45-
# Multiply constraint term with regular pdf using ROOT.RooProdPdf
46-
# Specify in fitTo() that internal constraints on parameter f should be
47-
# used
45+
# Multiply constraint term with regular pdf using ROOT.RooProdPdf Specify in
46+
# fitTo() that internal constraints on parameter f should be used
4847

4948
# Multiply constraint with pdf
5049
modelc = ROOT.RooProdPdf("modelc", "model with constraint", [model, fconstraint])
5150

5251
# Fit model (without use of constraint term)
53-
r1 = model.fitTo(d, Save=True)
52+
r1 = model.fitTo(d, Save=True, PrintLevel=-1)
5453

5554
# Fit modelc with constraint term on parameter f
56-
r2 = modelc.fitTo(d, Constrain={f}, Save=True)
55+
r2 = modelc.fitTo(d, Constrain={f}, Save=True, PrintLevel=-1)
5756

5857
# Method 2 - specify external constraint when fitting
5958
# ------------------------------------------------------------------------------------------
6059

6160
# Construct another Gaussian constraint pdf on parameter f at 0.8 with
6261
# resolution of 0.1
63-
fconstext = ROOT.RooGaussian("fconstext", "fconstext", f, ROOT.RooFit.RooConst(0.2), ROOT.RooFit.RooConst(0.1))
62+
fconstext = ROOT.RooGaussian("fconstext", "fconstext", ROOT.RooFit.RooConst(0.2), f, ROOT.RooFit.RooConst(0.1))
6463

6564
# Fit with external constraint
66-
r3 = model.fitTo(d, ExternalConstraints={fconstext}, Save=True)
65+
r3 = model.fitTo(d, ExternalConstraints={fconstext}, Save=True, PrintLevel=-1)
6766

6867
# Print the fit results
6968
print("fit result without constraint (data generated at f=0.5)")

0 commit comments

Comments
 (0)