Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooLagrangianMorphFunc.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * authors: *
5 * Lydia Brenner (lbrenner@cern.ch), Carsten Burgard (cburgard@cern.ch) *
6 * Katharina Ecker (kecker@cern.ch), Adam Kaluza (akaluza@cern.ch) *
7 * Copyright (c) 2000-2007, Regents of the University of California *
8 * and Stanford University. All rights reserved. *
9 * *
10 * Redistribution and use in source and binary forms, *
11 * with or without modification, are permitted according to the terms *
12 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
13 *****************************************************************************/
14
15////////////////////////////////////////////////////////////////////////////////////////////////
16//
17// RooLagrangianMorphFunc
18//
19// The RooLagrangianMorphFunc is a type of RooAbsReal that allows to morph
20// different input EFT samples to some arbitrary output EFT
21// sample, as long as the desired set of output parameters lie
22// within the realm spanned by the input samples. More
23// specifically, it expects as an input a TFile (or TDirectory)
24// with the following layout:
25//
26// TDirectory
27// |-sample1
28// | |-param_card // TH1 EFT parameter values of sample1
29// | | histogram1 // TH1 of some physics distribution
30// | |-subfolder1 // a subfolder (optional)
31// | | |-histogram2 // TH1 of some physics distribution
32// | | |-....
33// |-sample2
34// | |-param_card // TH1 of EFT parameter values of sample2
35// | | histogram1 // TH1 of some physics distribution
36// | |-subfolder1 // same folder structure as before
37// | | |-histogram2 // TH1 of some physics distribution
38// | | |-....
39// |-sampleN
40// The RooLagrangianMorphFunc operates on this structure, extracts data
41// and meta-data and produces a morphing result as a RooRealSumFunc
42// consisting of the input histograms with appropriate prefactors.
43//
44// The histograms to be morphed can be accessed via their paths in
45// the respective sample, e.g. using
46// "histogram"
47// or "subfolder1/histogram1"
48// or "some/deep/path/to/some/subfolder/histname"
49//
50////////////////////////////////////////////////////////////////////////////////////////////////
51
52#ifndef ROO_LAGRANGIAN_MORPH
53#define ROO_LAGRANGIAN_MORPH
54
55#include "RooFit/Floats.h"
56#include "RooAbsArg.h"
57#include "RooAbsReal.h"
58#include "RooArgList.h"
59#include "RooRatio.h"
60#include "RooRealSumFunc.h"
61#include "RooRealSumPdf.h"
62#include "RooSetProxy.h"
63#include "RooWrapperPdf.h"
64#include "TMatrixD.h"
65
66class RooWorkspace;
68class RooProduct;
69class RooRealVar;
70class TPair;
71class TFolder;
73
74#include <fstream>
75#include <iostream>
76#include <string>
77#include <vector>
78
80
81public:
82 typedef std::map<const std::string, double> ParamSet;
83 typedef std::map<const std::string, int> FlagSet;
84 typedef std::map<const std::string, ParamSet> ParamMap;
85 typedef std::map<const std::string, FlagSet> FlagMap;
86
87 struct Config {
88
89 std::string observableName;
90 std::string fileName;
93 std::vector<std::string> folderNames;
98 std::vector<RooArgList *> vertices;
99 std::vector<std::vector<const char *>> nonInterfering;
101 };
102
104 RooLagrangianMorphFunc(const char *name, const char *title, const char *filename, const char *observableName,
105 const RooArgSet &couplings, const RooArgSet &inputs);
106 RooLagrangianMorphFunc(const char *name, const char *title, const Config &config);
107 RooLagrangianMorphFunc(const RooLagrangianMorphFunc &other, const char *newName);
108
109 ~RooLagrangianMorphFunc() override;
110
111 std::list<double> *binBoundaries(RooAbsRealLValue & /*obs*/, double /*xlo*/, double /*xhi*/) const override;
112 std::list<double> *plotSamplingHint(RooAbsRealLValue & /*obs*/, double /*xlo*/, double /*xhi*/) const override;
113 bool isBinnedDistribution(const RooArgSet &obs) const override;
114 double evaluate() const override;
115 TObject *clone(const char *newname) const override;
116
117 bool checkObservables(const RooArgSet *nset) const override;
118 bool forceAnalyticalInt(const RooAbsArg &arg) const override;
119 Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet,
120 const char *rangeName = 0) const override;
121 double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName = 0) const override;
122 void printMetaArgs(std::ostream &os) const override;
123 RooAbsArg::CacheMode canNodeBeCached() const override;
124 void setCacheAndTrackHints(RooArgSet &) override;
125
127
128 void setParameters(const char *foldername);
129 void setParameters(TH1 *paramhist);
130 void setParameter(const char *name, double value);
131 void setFlag(const char *name, double value);
132 void setParameters(const ParamSet &params);
133 void setParameters(const RooArgList *list);
134 double getParameterValue(const char *name) const;
135 RooRealVar *getParameter(const char *name) const;
136 RooRealVar *getFlag(const char *name) const;
137 bool hasParameter(const char *paramname) const;
138 bool isParameterUsed(const char *paramname) const;
139 bool isParameterConstant(const char *paramname) const;
140 void setParameterConstant(const char *paramname, bool constant) const;
141 void setParameter(const char *name, double value, double min, double max);
142 void setParameter(const char *name, double value, double min, double max, double error);
143 void randomizeParameters(double z);
144 const RooArgSet *getParameterSet() const;
145 ParamSet getMorphParameters(const char *foldername) const;
147
149
150 int nParameters() const;
151 int nPolynomials() const;
152
153 bool isCouplingUsed(const char *couplname);
154 const RooArgList *getCouplingSet() const;
155 ParamSet getCouplings() const;
156
157 TMatrixD getMatrix() const;
159 double getCondition() const;
160
161 RooRealVar *getObservable() const;
162 RooRealVar *getBinWidth() const;
163
164 void printEvaluation() const;
165 void printCouplings() const;
166 void printFlags() const;
167 void printPhysics() const;
168
169 RooProduct *getSumElement(const char *name) const;
170
171 std::vector<std::string> getSamples() const;
172
173 double expectedUncertainty() const;
174 TH1 *createTH1(const std::string &name);
175 TH1 *createTH1(const std::string &name, bool correlateErrors);
176
177private:
178 class CacheElem;
179 void init();
180 void setup(bool ownParams = true);
181 void disableInterference(const std::vector<const char *> &nonInterfering);
182 void disableInterferences(const std::vector<std::vector<const char *>> &nonInterfering);
183
184 void addFolders(const RooArgList &folders);
185
186 bool hasCache() const;
188 void updateSampleWeights();
189
190 RooRealVar *setupObservable(const char *obsname, TClass *mode, TObject *inputExample);
191
192public:
193 /// length of floating point digits precision supported by implementation.
194 static constexpr double implementedPrecision = RooFit::SuperFloatPrecision::digits10;
195
196 void writeMatrixToFile(const TMatrixD &matrix, const char *fname);
197 void writeMatrixToStream(const TMatrixD &matrix, std::ostream &stream);
198 TMatrixD readMatrixFromFile(const char *fname);
199 TMatrixD readMatrixFromStream(std::istream &stream);
200
201 int countSamples(std::vector<RooArgList *> &vertices);
202 int countSamples(int nprod, int ndec, int nboth);
203
204 std::map<std::string, std::string>
205 createWeightStrings(const ParamMap &inputs, const std::vector<std::vector<std::string>> &vertices);
206 std::map<std::string, std::string>
207 createWeightStrings(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings);
208 std::map<std::string, std::string>
209 createWeightStrings(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings,
210 const FlagMap &flagValues, const RooArgList &flags,
211 const std::vector<RooArgList *> &nonInterfering);
212 RooArgSet createWeights(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings,
213 const FlagMap &inputFlags, const RooArgList &flags,
214 const std::vector<RooArgList *> &nonInterfering);
215 RooArgSet createWeights(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings);
216
217 bool updateCoefficients();
218 bool useCoefficients(const TMatrixD &inverse);
219 bool useCoefficients(const char *filename);
220 bool writeCoefficients(const char *filename);
221
222 int countContributingFormulas() const;
223 RooAbsReal *getSampleWeight(const char *name);
224 void printParameters(const char *samplename) const;
225 void printParameters() const;
226 void printSamples() const;
227 void printSampleWeights() const;
228 void printWeights() const;
229
230 void setScale(double val);
231 double getScale();
232
233 int nSamples() const { return _config.folderNames.size(); }
234
235 RooRealSumFunc *getFunc() const;
236 std::unique_ptr<RooWrapperPdf> createPdf() const;
237
239 double expectedEvents(const RooArgSet *nset) const;
240 double expectedEvents(const RooArgSet &nset) const;
241 double expectedEvents() const;
242 bool selfNormalized() const { return true; }
243
246
247 static std::unique_ptr<RooRatio> makeRatio(const char *name, const char *title, RooArgList &nr, RooArgList &dr);
248
249private:
250 mutable RooObjCacheManager _cacheMgr; //! The cache manager
251 double _scale = 1.0;
252 std::map<std::string, int> _sampleMap;
259 std::vector<std::vector<RooListProxy *>> _diagrams;
260 std::vector<RooListProxy *> _nonInterfering;
261
263};
264
265#endif
#define f(i)
Definition RSha256.hxx:104
#define ClassDefOverride(name, id)
Definition Rtypes.h:329
char name[80]
Definition TGX11.cxx:110
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:35
Class RooLagrangianMorphing is a implementation of the method of Effective Lagrangian Morphing,...
static constexpr double implementedPrecision
length of floating point digits precision supported by implementation.
bool isParameterConstant(const char *paramname) const
return true if the parameter with the given name is set constant, false otherwise
bool isBinnedDistribution(const RooArgSet &obs) const override
check if this PDF is a binned distribution in the given observable
int nPolynomials() const
return the number of samples in this morphing function
void setParameter(const char *name, double value)
set one parameter to a specific value
ParamSet getMorphParameters() const
retrieve the parameter set
double evaluate() const override
call getVal on the internal function
void disableInterference(const std::vector< const char * > &nonInterfering)
disable interference between terms
RooProduct * getSumElement(const char *name) const
return the RooProduct that is the element of the RooRealSumPdfi corresponding to the given sample nam...
void insert(RooWorkspace *ws)
RooRealVar * getBinWidth() const
retrieve the histogram observable
void writeMatrixToFile(const TMatrixD &matrix, const char *fname)
write a matrix to a file
RooRealVar * getParameter(const char *name) const
retrieve the RooRealVar object incorporating the parameter with the given name
bool useCoefficients(const TMatrixD &inverse)
setup the morphing function with a predefined inverse matrix call this function before any other afte...
const RooArgSet * getParameterSet() const
get the set of parameters
std::map< const std::string, int > FlagSet
TMatrixD readMatrixFromStream(std::istream &stream)
read a matrix from a stream
std::vector< std::string > getSamples() const
return the vector of sample names, used to build the morph func
int countSamples(std::vector< RooArgList * > &vertices)
calculate the number of samples needed to morph a certain physics process
void setCacheAndTrackHints(RooArgSet &) override
Retrieve the matrix of coefficients.
ParamSet getCouplings() const
retrieve a set of couplings (-?-)
void printSampleWeights() const
print the current sample weights
std::map< const std::string, double > ParamSet
RooArgSet createWeights(const ParamMap &inputs, const std::vector< RooArgList * > &vertices, RooArgList &couplings, const FlagMap &inputFlags, const RooArgList &flags, const std::vector< RooArgList * > &nonInterfering)
create only the weight formulas. static function for external usage.
RooLagrangianMorphFunc * getLinear() const
void writeMatrixToStream(const TMatrixD &matrix, std::ostream &stream)
write a matrix to a stream
void addFolders(const RooArgList &folders)
convert the RooArgList folders into a simple vector of std::string
std::map< const std::string, ParamSet > ParamMap
bool updateCoefficients()
retrive the new physics objects and update the weights in the morphing function
double _scale
The cache manager.
RooRealVar * getObservable() const
retrieve the histogram observable
int countContributingFormulas() const
count the number of formulas that correspond to the current parameter set
std::list< double > * plotSamplingHint(RooAbsRealLValue &, double, double) const override
retrieve the sample Hint
RooRealVar * getFlag(const char *name) const
retrieve the RooRealVar object incorporating the flag with the given name
void randomizeParameters(double z)
randomize the parameters a bit useful to test and debug fitting
bool isCouplingUsed(const char *couplname)
check if there is any morphing power provided for the given coupling morphing power is provided as so...
void readParameters(TDirectory *f)
read the parameters from the input file
double getScale()
get energy scale of the EFT expansion
std::vector< RooListProxy * > _nonInterfering
double getCondition() const
Retrieve the condition of the coefficient matrix.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const override
Retrieve the matrix of coefficients.
TMatrixD getMatrix() const
Retrieve the matrix of coefficients.
void printWeights() const
print the current sample weights
void printCouplings() const
print a set of couplings
TMatrixD readMatrixFromFile(const char *fname)
read a matrix from a text file
~RooLagrangianMorphFunc() override
default destructor
void printParameters() const
print the parameters and their current values
void printPhysics() const
print the current physics values
static std::unique_ptr< RooRatio > makeRatio(const char *name, const char *title, RooArgList &nr, RooArgList &dr)
Return the RooRatio form of products and denominators of morphing functions.
void setFlag(const char *name, double value)
set one flag to a specific value
TH1 * createTH1(const std::string &name)
retrieve a histogram output of the current morphing settings
double expectedUncertainty() const
return the expected uncertainty for the current parameter set
int nParameters() const
return the number of parameters in this morphing function
bool hasParameter(const char *paramname) const
check if a parameter of the given name is contained in the list of known parameters
bool checkObservables(const RooArgSet *nset) const override
check if observable exists in the RooArgSet (-?-)
bool hasCache() const
return true if a cache object is present, false otherwise
void printFlags() const
print the flags and their current values
void setScale(double val)
set energy scale of the EFT expansion
TMatrixD getInvertedMatrix() const
Retrieve the matrix of coefficients after inversion.
RooAbsArg::CacheMode canNodeBeCached() const override
Retrieve the matrix of coefficients.
void updateSampleWeights()
update sample weight (-?-)
void setParameters(const char *foldername)
set the morphing parameters to those supplied in the sample with the given name
bool isParameterUsed(const char *paramname) const
check if there is any morphing power provided for the given parameter morphing power is provided as s...
RooAbsPdf::ExtendMode extendMode() const
return extended mored capabilities
std::map< const std::string, FlagSet > FlagMap
bool forceAnalyticalInt(const RooAbsArg &arg) const override
Force analytical integration for the given observable.
void setParameterConstant(const char *paramname, bool constant) const
call setConstant with the boolean argument provided on the parameter with the given name
void disableInterferences(const std::vector< std::vector< const char * > > &nonInterfering)
disable interference between terms
void printSamples() const
print all the known samples to the console
double getParameterValue(const char *name) const
set one parameter to a specific value
void setup(bool ownParams=true)
setup this instance with the given set of operators and vertices if own=true, the class will own the ...
void printMetaArgs(std::ostream &os) const override
Retrieve the matrix of coefficients.
std::unique_ptr< RooWrapperPdf > createPdf() const
(currently similar to cloning the Pdf
RooAbsReal * getSampleWeight(const char *name)
retrieve the weight (prefactor) of a sample with the given name
std::map< std::string, std::string > createWeightStrings(const ParamMap &inputs, const std::vector< std::vector< std::string > > &vertices)
create only the weight formulas. static function for external usage.
std::vector< std::vector< RooListProxy * > > _diagrams
double expectedEvents() const
return the number of expected events for the current parameter set
std::map< std::string, int > _sampleMap
TObject * clone(const char *newname) const override
cloning method
RooLagrangianMorphFunc::CacheElem * getCache() const
retrieve the cache object
RooRealVar * setupObservable(const char *obsname, TClass *mode, TObject *inputExample)
setup observable, recycle existing observable if defined
const RooArgList * getCouplingSet() const
get the set of couplings
RooRealSumFunc * getFunc() const
get the func
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
retrieve the list of bin boundaries
void printEvaluation() const
print the contributing samples and their respective weights
bool writeCoefficients(const char *filename)
write the inverse matrix to a file
void collectInputs(TDirectory *f)
retrieve the physics inputs
void init()
initialise inputs required for the morphing function
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=0) const override
Retrieve the mat.
RooListProxy is the concrete proxy for RooArgList objects.
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
A histogram function that assigns scale parameters to every bin.
A RooProduct represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
RooSetProxy is the concrete proxy for RooArgSet objects.
Definition RooSetProxy.h:23
The RooWorkspace is a persistable container for RooFit projects.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:80
Describe directory structure in memory.
Definition TDirectory.h:45
A TFolder object is a collection of objects and folders.
Definition TFolder.h:30
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:58
Mother of all ROOT objects.
Definition TObject.h:41
Class used by TMap to store (key,value) pairs.
Definition TMap.h:102
std::vector< std::string > folderNames
std::vector< std::vector< const char * > > nonInterfering
std::vector< RooArgList * > vertices
void ws()
Definition ws.C:66