1#ifndef INCLUDE_HISTFACTORYNAVIGATION_H
2#define INCLUDE_HISTFACTORYNAVIGATION_H
27 const std::string& WorkspaceName=
"combined",
28 const std::string& ModelConfigName=
"ModelConfig");
43 bool IncludeConstantParams=
false);
47 bool IncludeConstantParams=
false);
60 double GetBinValue(
int bin,
const std::string& channel);
62 double GetBinValue(
int bin,
const std::string& channel,
const std::string& sample);
67 const std::string& sample,
const std::string&
name=
"");
100 void SetConstant(
const std::string& regExpr=
".*",
bool constant=
true);
#define ClassDef(name, id)
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Abstract interface for all probability density functions.
Abstract base class for objects that represent a real value and implements functionality common to al...
RooArgList is a container object that can hold multiple RooAbsArg objects.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
RooDataSet is a container class to hold unbinned data.
A RooProduct represents the product of a given set of RooAbsReal objects.
RooRealVar represents a variable that can be changed from the outside.
void PrintParameters(bool IncludeConstantParams=false)
Print the current values and errors of pdf parameters.
std::vector< std::string > fChannelNameVec
The list of channels.
void SetMinBinToPrint(int min)
RooArgSet * GetObservableSet(const std::string &channel)
Get the set of observables for a given channel.
int GetMinBinToPrint() const
TH1 * GetDataHist(RooDataSet *data, const std::string &channel, const std::string &name="")
Get a histogram from the dataset for this channel.
RooArgSet * fObservables
The observables.
void _GetNodes(ModelConfig *mc)
Fetch the node information for the pdf in question, and save it in the various collections in this cl...
std::map< std::string, RooAbsPdf * > fChannelSumNodeMap
Map of channel names to pdf without constraint.
void ReplaceNode(const std::string &ToReplace, RooAbsArg *ReplaceWith)
Find a node in the pdf and replace it with a new node These nodes can be functions,...
void PrintState()
Should pretty print all channels and the current values
void PrintSampleComponents(const std::string &channel, const std::string &sample)
Print the different components that make up a sample (NormFactors, Statistical Uncertainties,...
void SetPrintWidths(const std::string &channel)
Set the title and bin widths.
RooAbsPdf * fModel
The HistFactory Pdf Pointer.
TH1 * GetChannelHist(const std::string &channel, const std::string &name="")
Get the total channel histogram for this channel.
TH1 * GetSampleHist(const std::string &channel, const std::string &sample, const std::string &name="")
The (current) histogram for that sample This includes all parameters and interpolation.
TH1 * MakeHistFromRooFunction(RooAbsReal *func, RooArgList vars, std::string name="Hist")
Make a histogram from a function Edit so it can take a RooArgSet of parameters.
void PrintMultiDimHist(TH1 *hist, int bin_print_width)
Print a histogram's contents to the screen void PrettyPrintHistogram(TH1* hist);.
virtual ~HistFactoryNavigation()
void PrintSampleParameters(const std::string &channel, const std::string &sample, bool IncludeConstantParams=false)
Print parameters that effect a particular sample.
std::map< std::string, RooAbsReal * > GetSampleFunctionMap(const std::string &channel)
Get a map of sample names to their functions for a particular channel.
RooAbsReal * SampleFunction(const std::string &channel, const std::string &sample)
Get the RooAbsReal function for a given sample in a given channel.
double GetBinValue(int bin, const std::string &channel)
The value of the ith bin for the total in that channel.
RooAbsArg * findChild(const std::string &name, RooAbsReal *parent) const
Internal method implementation of finding a daughter node from a parent node (looping over all genera...
std::map< std::string, std::map< std::string, RooAbsReal * > > fChannelSampleFunctionMap
Map of Map of Channel, Sample names to Function Nodes Used by doing: fChannelSampleFunctionMap["MyCha...
void SetMaxBinToPrint(int max)
std::map< std::string, RooAbsPdf * > fChannelPdfMap
Map of channel names to their full pdf's.
std::vector< std::string > GetChannelSampleList(const std::string &channel)
RooAbsPdf * GetModel() const
Get the model for this channel.
void PrintModelAndData(RooDataSet *data)
Print the model and the data, comparing channel by channel.
std::map< std::string, RooArgSet * > fChannelObservMap
Map of channel names to their set of ovservables.
RooAbsPdf * GetChannelPdf(const std::string &channel)
void PrintDataSet(RooDataSet *data, const std::string &channel="")
Print a "HistFactory style" RooDataSet in a readable way.
void PrintChannelParameters(const std::string &channel, bool IncludeConstantParams=false)
Print parameters that effect a particular channel.
void SetConstant(const std::string ®Expr=".*", bool constant=true)
RooArgSet _GetAllProducts(RooProduct *node)
Recursively get all products of products.
THStack * GetChannelStack(const std::string &channel, const std::string &name="")
Get a stack of all samples in a channel.
int GetMaxBinToPrint() const
RooAbsReal * GetConstraintTerm(const std::string ¶meter)
Get the constraint term for a given systematic (alpha or gamma)
void DrawChannel(const std::string &channel, RooDataSet *data=nullptr)
Draw a stack of the channel, and include data if the pointer is supplied.
RooRealVar * var(const std::string &varName) const
double GetConstraintUncertainty(const std::string ¶meter)
Get the uncertainty based on the constraint term for a given systematic.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
TH1 is the base class of all histogram classes in ROOT.
The Histogram stack class.
RooStats::ModelConfig ModelConfig
Namespace for the RooStats classes.