Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ModifyInterpolation.C
Go to the documentation of this file.
1/*
2An example script to modify the interpolation used in HistFactory models.
3Usage:
4
5make the standard example workspace with histfactory
6$ prepareHistFactory
7$ hist2workspace config/example.xml
8
9inspect file (note we are using the model with Gaussian constraints)
10$ root.exe results/example_combined_GaussExample_model.root
11root [1] combined->Print()
12
13notice there is a new class in the new version:
14RooStats::HistFactory::FlexibleInterpVar::background1_channel1_epsilon[ paramList=(alpha_syst2) ] = 0.999209
15
16You only need this macro if you don't want the CMS style normalization interpolation.
17
18root [2] .L $ROOTSYS/tutorials/histfactory/ModifyInterpolation.C+
19root [3] CheckInterpolation(combined)
20[#1] INFO:InputArguments -- interp code for alpha_syst1 = 1
21[#1] INFO:InputArguments -- interp code for alpha_syst2 = 1
22[#1] INFO:InputArguments -- interp code for alpha_syst3 = 1
23
24One can change the interpolation for a specific set of parameters
25root [5] combined->set("ModelConfig_NuisParams")->Print()
26RooArgSet:: = (alpha_syst2,alpha_syst3)
27root [6] ModifyInterpolationForSet(combined->set("ModelConfig_NuisParams"),2)
28root [7] CheckInterpolation(combined)
29[#1] INFO:InputArguments -- interp code for alpha_syst1 = 1
30[#1] INFO:InputArguments -- interp code for alpha_syst2 = 2
31[#1] INFO:InputArguments -- interp code for alpha_syst3 = 2
32
33Or one can change the interpolation type for all the parameters:
34root [11] ModifyInterpolationForAll(combined,2)
35root [12] CheckInterpolation(combined)
36[#1] INFO:InputArguments -- interp code for alpha_syst1 = 2
37[#1] INFO:InputArguments -- interp code for alpha_syst2 = 2
38[#1] INFO:InputArguments -- interp code for alpha_syst3 = 2
39
40and then you can save the workspace with those modifications
41root [13] combined->writeToFile("testModified.root")
42
43then test to make sure the changes are in the file
44$ root.exe testModified.root
45root [1] .L $ROOTSYS/tutorials/histfactory/ModifyInterpolation.C+
46root [2] CheckInterpolation(combined)
47[#1] INFO:InputArguments -- interp code for alpha_syst1 = 2
48[#1] INFO:InputArguments -- interp code for alpha_syst2 = 2
49[#1] INFO:InputArguments -- interp code for alpha_syst3 = 2
50*/
51#include "RooRealVar.h"
52#include "TIterator.h"
53#include "RooWorkspace.h"
55
56using namespace RooStats;
57using namespace HistFactory;
58
59// define functions for ACLIC
60void ModifyInterpolationForAll(RooWorkspace* ws, int code=1);
61void ModifyInterpolationForSet(RooArgSet* modifySet, int code = 1);
63
64// Codes for interpolation
65// code = 0: piece-wise linear
66// code = 1: pice-wise log
67// code = 2: parabolic interp with linear extrap
68// code = 3: parabolic version of log-normal
69
71 cout <<"Choose from the following"<<endl;
72 cout <<"void ModifyInterpolationForAll(RooWorkspace* ws, int code=1);"<<endl;
73 cout <<"void ModifyInterpolationForSet(RooArgSet* modifySet, int code = 1);"<<endl;
74 cout <<"void CheckInterpolation(RooWorkspace* ws);"<<endl;
75}
76
79 TIter it = funcs.createIterator();
80 TObject* tempObj = nullptr;
81 while((tempObj=it.Next())){
82 FlexibleInterpVar* flex = dynamic_cast<FlexibleInterpVar*>(tempObj);
83 if(flex){
84 flex->setAllInterpCodes(code);
85 }
86 }
87}
88
89void ModifyInterpolationForSet(RooArgSet* modifySet, int code){
90
91 TIter it = modifySet->createIterator();
92 RooRealVar* alpha = nullptr;
93 while((alpha=(RooRealVar*)it.Next())){
94 TIter serverIt = alpha->clientIterator();
95 TObject* tempObj = nullptr;
96 while((tempObj=serverIt.Next())){
97 FlexibleInterpVar* flex = dynamic_cast<FlexibleInterpVar*>(tempObj);
98 if(flex){
99 flex->printAllInterpCodes();
100 flex->setInterpCode(*alpha,code);
101 flex->printAllInterpCodes();
102 }
103 }
104 }
105
106}
107
108
111 TIter it = funcs.createIterator();
112 TObject* tempObj=0;
113 while((tempObj=it.Next())){
114 FlexibleInterpVar* flex = dynamic_cast<FlexibleInterpVar*>(tempObj);
115 if(flex){
116 flex->printAllInterpCodes();
117 }
118 }
119}
void CheckInterpolation(RooWorkspace *ws)
void ModifyInterpolation()
void ModifyInterpolationForSet(RooArgSet *modifySet, int code=1)
void ModifyInterpolationForAll(RooWorkspace *ws, int code=1)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void funcs
TIterator * clientIterator() const
Retrieve a client iterator.
Definition RooAbsArg.h:133
TIterator * createIterator(bool dir=kIterForward) const
TIterator-style iteration over contained elements.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setInterpCode(RooAbsReal &param, int code)
The RooWorkspace is a persistable container for RooFit projects.
RooArgSet allFunctions() const
Return set with all function objects.
TObject * Next()
Mother of all ROOT objects.
Definition TObject.h:41
Namespace for the RooStats classes.
Definition Asimov.h:19