Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
7 changes: 7 additions & 0 deletions .github/workflows/cvmfs-ci.yml
Original file line number Diff line number Diff line change
Expand Up @@ -283,6 +283,13 @@ jobs:
cd test
python3 test_interference.py

- uses: ./.github/actions/run-in-cvmfs
name: FastScan template analysis CMSHistFunc
with:
script: |
text2workspace.py data/ci/template-analysis_shapeInterp.txt -o ws_template-analysis.root --mass 200
combineTool.py -M FastScan -w ws_template-analysis.root:w

- name: Prepare code coverage report
if: ${{ success() && (matrix.coverage == true) }}
run: |
Expand Down
8 changes: 2 additions & 6 deletions interface/Combine.h
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@
#include "RooAbsReal.h"
#include "RooRealVar.h"

#include "CombineUtils.h"

class TDirectory;
class TTree;
class LimitAlgo;
Expand Down Expand Up @@ -123,10 +125,4 @@ class Combine {
static std::string textToWorkspaceString_;
};

std::unique_ptr<RooAbsReal> combineCreateNLL(RooAbsPdf &pdf,
RooAbsData &data,
RooArgSet const *constraint = nullptr,
bool offset = false,
bool warnAboutDifferentBackend = true);

#endif
29 changes: 29 additions & 0 deletions interface/CombineUtils.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,29 @@
#ifndef HiggsAnalysis_CombinedLimit_CombineUtils_h
#define HiggsAnalysis_CombinedLimit_CombineUtils_h

#include <memory>

class RooAbsReal;
class RooAbsPdf;
class RooAbsData;
class RooArgSet;

std::unique_ptr<RooAbsReal> combineCreateNLL(RooAbsPdf& pdf,
RooAbsData& data,
RooArgSet const* constraint = nullptr,
bool offset = false,
bool warnAboutDifferentBackend = true);

// Class that makes the utility functions also available in Python (adding classes to the dictionaries is less brittle than adding functions)
class CombineUtils {
public:
static std::unique_ptr<RooAbsReal> combineCreateNLL(RooAbsPdf& pdf,
RooAbsData& data,
RooArgSet const* constraint = nullptr,
bool offset = false,
bool warnAboutDifferentBackend = true) {
return ::combineCreateNLL(pdf, data, constraint, offset, warnAboutDifferentBackend);
}
};

#endif
7 changes: 4 additions & 3 deletions python/tool_base/FastScan.py
Original file line number Diff line number Diff line change
Expand Up @@ -57,7 +57,7 @@ def run_method(self):
data = f_d.Get(ws_d[1])
else:
data = f_d.Get(ws_d[1]).data(ws_d[2])
nll = ROOT.combineCreateNLL(pdf, data)
nll = ROOT.CombineUtils.combineCreateNLL(pdf, data)
pars = pdf.getParameters(data)
pars.Print()
snap = pars.snapshot()
Expand Down Expand Up @@ -114,7 +114,8 @@ def run_method(self):
grd1.Write()
grd2.Write()
pars.assignValueOnly(snap)
canv = ROOT.TCanvas(self.args.output, self.args.output)
canv_name = "par_%s" % par.GetName()
canv = ROOT.TCanvas(canv_name, canv_name)
pads = plot.MultiRatioSplit([0.4, 0.3], [0.005, 0.005], [0.005, 0.005])
pads[0].cd()
plot.Set(gr, MarkerSize=0.5)
Expand Down Expand Up @@ -143,7 +144,7 @@ def run_method(self):
if page == len(doPars) - 1:
extra = ")"
print(extra)
canv.Print(".pdf%s" % extra)
canv.Print("%s.pdf%s" % (self.args.output, extra))
page += 1

outfile.Write()
2 changes: 1 addition & 1 deletion python/tool_base/TaylorExpand.py
Original file line number Diff line number Diff line change
Expand Up @@ -186,7 +186,7 @@ def load_workspace(self, file, POIs, data="data_obs", snapshot="MultiDimFit"):
mc = self.loaded_wsp.genobj("ModelConfig")
pdf = mc.GetPdf()
data = self.loaded_wsp.data(data)
self.nll = ROOT.combineCreateNLL(pdf, data)
self.nll = ROOT.CombineUtils.combineCreateNLL(pdf, data)
self.loaded_wsp.loadSnapshot("MultiDimFit")
print("...NLL loaded")
# nll.Print()
Expand Down
1 change: 1 addition & 0 deletions src/classes.h
Original file line number Diff line number Diff line change
Expand Up @@ -63,3 +63,4 @@
#include "HiggsAnalysis/CombinedLimit/interface/RooEFTScalingFunction.h"

#include "HiggsAnalysis/CombinedLimit/interface/CombineCodegenImpl.h"
#include "HiggsAnalysis/CombinedLimit/interface/CombineUtils.h"
1 change: 1 addition & 0 deletions src/classes_def.xml
Original file line number Diff line number Diff line change
Expand Up @@ -220,4 +220,5 @@
<class name="RooSumTwoExpPdf" />
<class name="RooPowerLawPdf" />
<class name="RooSumTwoPowerLawPdf" />
<class name="CombineUtils" />
</lcgdict>
Loading