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 "RooAbsReal.h"
57#include "RooArgList.h"
58#include "RooRatio.h"
59#include "RooRealSumFunc.h"
60#include "RooRealSumPdf.h"
61#include "RooSetProxy.h"
62#include "RooWrapperPdf.h"
63#include "TMatrixD.h"
64
65class RooWorkspace;
66class RooProduct;
67class RooRealVar;
68class TPair;
69class TFolder;
71
72#include <fstream>
73#include <iostream>
74#include <string>
75#include <vector>
76
78
79public:
80 typedef std::map<const std::string, double> ParamSet;
81 typedef std::map<const std::string, int> FlagSet;
82 typedef std::map<const std::string, ParamSet> ParamMap;
83 typedef std::map<const std::string, FlagSet> FlagMap;
84
85 struct Config {
86 std::string observableName;
88 std::string fileName;
91 std::vector<std::string> folderNames;
96 std::vector<RooArgList *> vertices;
97 std::vector<std::vector<const char *>> nonInterfering;
99 bool normalize = false;
100 };
101
103 RooLagrangianMorphFunc(const char *name, const char *title, const Config &config);
104 RooLagrangianMorphFunc(const RooLagrangianMorphFunc &other, const char *newName);
105
106 ~RooLagrangianMorphFunc() override;
107
108 std::list<double> *binBoundaries(RooAbsRealLValue & /*obs*/, double /*xlo*/, double /*xhi*/) const override;
109 std::list<double> *plotSamplingHint(RooAbsRealLValue & /*obs*/, double /*xlo*/, double /*xhi*/) const override;
110 bool isBinnedDistribution(const RooArgSet &obs) const override;
111 double evaluate() const override;
112 TObject *clone(const char *newname) const override { return new RooLagrangianMorphFunc(*this, newname); }
113
114 bool checkObservables(const RooArgSet *nset) const override;
115 bool forceAnalyticalInt(const RooAbsArg &arg) const override;
116 Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet,
117 const char *rangeName = nullptr) const override;
118 double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName = nullptr) const override;
119 void printMetaArgs(std::ostream &os) const override;
120 RooAbsArg::CacheMode canNodeBeCached() const override;
121 void setCacheAndTrackHints(RooArgSet &) override;
122
124
125 void setParameters(const char *foldername);
126 void setParameters(TH1 *paramhist);
127 void setParameter(const char *name, double value);
128 void setFlag(const char *name, double value);
129 void setParameters(const ParamSet &params);
130 void setParameters(const RooArgList *list);
131 double getParameterValue(const char *name) const;
132 RooRealVar *getParameter(const char *name) const;
133 RooRealVar *getFlag(const char *name) const;
134 bool hasParameter(const char *paramname) const;
135 bool isParameterUsed(const char *paramname) const;
136 bool isParameterConstant(const char *paramname) const;
137 void setParameterConstant(const char *paramname, bool constant) const;
138 void setParameter(const char *name, double value, double min, double max);
139 void setParameter(const char *name, double value, double min, double max, double error);
140 void randomizeParameters(double z);
141 const RooArgSet *getParameterSet() const;
142 ParamSet getMorphParameters(const char *foldername) const;
144
146
147 int nParameters() const;
148 int nPolynomials() const;
149
150 bool isCouplingUsed(const char *couplname);
151 const RooArgList *getCouplingSet() const;
152 ParamSet getCouplings() const;
153
154 TMatrixD getMatrix() const;
156 double getCondition() const;
157
158 RooRealVar *getObservable() const;
159 RooRealVar *getBinWidth() const;
160
161 void printEvaluation() const;
162 void printCouplings() const;
163 void printFlags() const;
164 void printPhysics() const;
165
166 RooProduct *getSumElement(const char *name) const;
167
168 std::vector<std::string> getSamples() const;
169
170 double expectedUncertainty() const;
171 TH1 *createTH1(const std::string &name);
172 TH1 *createTH1(const std::string &name, bool correlateErrors);
173
174private:
175 class CacheElem;
176 void init();
177 void setup(bool ownParams = true);
178 void disableInterference(const std::vector<const char *> &nonInterfering);
179 void disableInterferences(const std::vector<std::vector<const char *>> &nonInterfering);
180
181 bool hasCache() const;
183 void updateSampleWeights();
184
185 RooRealVar *setupObservable(const char *obsname, TClass *mode, TObject *inputExample);
186
187public:
188 /// length of floating point digits precision supported by implementation.
189 static constexpr double implementedPrecision = RooFit::SuperFloatPrecision::digits10;
190
191 void writeMatrixToFile(const TMatrixD &matrix, const char *fname);
192 void writeMatrixToStream(const TMatrixD &matrix, std::ostream &stream);
193 TMatrixD readMatrixFromFile(const char *fname);
194 TMatrixD readMatrixFromStream(std::istream &stream);
195
196 int countSamples(std::vector<RooArgList *> &vertices);
197 int countSamples(int nprod, int ndec, int nboth);
198
199 std::map<std::string, std::string>
200 createWeightStrings(const ParamMap &inputs, const std::vector<std::vector<std::string>> &vertices);
201 std::map<std::string, std::string>
202 createWeightStrings(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings);
203 std::map<std::string, std::string>
204 createWeightStrings(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings,
205 const FlagMap &flagValues, const RooArgList &flags,
206 const std::vector<std::vector<std::string>> &nonInterfering);
207 RooArgSet createWeights(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings,
208 const FlagMap &inputFlags, const RooArgList &flags,
209 const std::vector<std::vector<std::string>> &nonInterfering);
210 RooArgSet createWeights(const ParamMap &inputs, const std::vector<RooArgList *> &vertices, RooArgList &couplings);
211
212 bool updateCoefficients();
213 bool useCoefficients(const TMatrixD &inverse);
214 bool useCoefficients(const char *filename);
215 bool writeCoefficients(const char *filename);
216
217 int countContributingFormulas() const;
218 RooAbsReal *getSampleWeight(const char *name);
219 void printParameters(const char *samplename) const;
220 void printParameters() const;
221 void printSamples() const;
222 void printSampleWeights() const;
223 void printWeights() const;
224
225 void setScale(double val);
226 double getScale();
227
228 int nSamples() const { return _config.folderNames.size(); }
229
230 RooRealSumFunc *getFunc() const;
231 std::unique_ptr<RooWrapperPdf> createPdf() const;
232
234 double expectedEvents(const RooArgSet *nset) const;
235 double expectedEvents(const RooArgSet &nset) const;
236 double expectedEvents() const;
237 bool selfNormalized() const { return true; }
238
241
242 static std::unique_ptr<RooRatio> makeRatio(const char *name, const char *title, RooArgList &nr, RooArgList &dr);
243
244private:
245 mutable RooObjCacheManager _cacheMgr; //! The cache manager
246 double _scale = 1.0;
247 std::map<std::string, int> _sampleMap;
254 std::vector<std::vector<RooListProxy *>> _diagrams;
255 std::vector<std::vector<std::string>> _nonInterfering;
256
258};
259
260#endif
#define f(i)
Definition RSha256.hxx:104
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void value
Option_t Option_t TPoint TPoint const char mode
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:74
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:62
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:55
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
RooArgSet createWeights(const ParamMap &inputs, const std::vector< RooArgList * > &vertices, RooArgList &couplings, const FlagMap &inputFlags, const RooArgList &flags, const std::vector< std::vector< std::string > > &nonInterfering)
create only the weight formulas. static function for external usage.
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::vector< std::string > > _nonInterfering
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
RooLagrangianMorphFunc * getLinear() const
void writeMatrixToStream(const TMatrixD &matrix, std::ostream &stream)
write a matrix to a stream
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
double getCondition() const
Retrieve the condition of the coefficient matrix.
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.
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Retrieve the matrix of coefficients.
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.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &numVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Retrieve the mat.
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
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
Class RooObjCacheManager is an implementation of class RooCacheManager<RooAbsCacheElement> and specia...
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:40
The RooWorkspace is a persistable container for RooFit projects.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
Describe directory structure in memory.
Definition TDirectory.h:45
<div class="legacybox"><h2>Legacy Code</h2> TFolder is a legacy interface: there will be no bug fixes...
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