ROOT logo
ROOT » ROOFIT » HISTFACTORY » RooStats::HistFactory::HistFactoryNavigation

class RooStats::HistFactory::HistFactoryNavigation

Function Members (Methods)

public:
virtual~HistFactoryNavigation()
static TClass*Class()
voidDrawChannel(const string& channel, RooDataSet* data = NULL)
doubleGetBinValue(int bin, const string& channel)
doubleGetBinValue(int bin, const string& channel, const string& sample)
TH1*GetChannelHist(const string& channel, const string& name = "")
RooAbsPdf*GetChannelPdf(const string& channel)
vector<std::string>GetChannelSampleList(const string& channel)
THStack*GetChannelStack(const string& channel, const string& name = "")
RooAbsReal*GetConstraintTerm(const string& parameter)
doubleGetConstraintUncertainty(const string& parameter)
TH1*GetDataHist(RooDataSet* data, const string& channel, const string& name = "")
intGetMaxBinToPrint() const
intGetMinBinToPrint() const
RooAbsPdf*GetModel() const
RooArgSet*GetObservableSet(const string& channel)
TH1*GetSampleHist(const string& channel, const string& sample, const string& name = "")
RooStats::HistFactory::HistFactoryNavigationHistFactoryNavigation(RooStats::ModelConfig* mc)
RooStats::HistFactory::HistFactoryNavigationHistFactoryNavigation(const RooStats::HistFactory::HistFactoryNavigation&)
RooStats::HistFactory::HistFactoryNavigationHistFactoryNavigation(RooAbsPdf* model, RooArgSet* observables)
RooStats::HistFactory::HistFactoryNavigationHistFactoryNavigation(const string& File, const string& WorkspaceName = "combined", const string& ModelConfigName = "ModelConfig")
virtual TClass*IsA() const
RooStats::HistFactory::HistFactoryNavigation&operator=(const RooStats::HistFactory::HistFactoryNavigation&)
voidPrintChannelParameters(const string& channel, bool IncludeConstantParams = false)
voidPrintDataSet(RooDataSet* data, const string& channel = "")
voidPrintModelAndData(RooDataSet* data)
voidPrintParameters(bool IncludeConstantParams = false)
voidPrintSampleComponents(const string& channel, const string& sample)
voidPrintSampleParameters(const string& channel, const string& sample, bool IncludeConstantParams = false)
voidPrintState()
voidPrintState(const string& channel)
voidReplaceNode(const string& ToReplace, RooAbsArg* ReplaceWith)
RooAbsReal*SampleFunction(const string& channel, const string& sample)
voidSetConstant(const string& regExpr = ".*", bool constant = true)
voidSetMaxBinToPrint(int max)
voidSetMinBinToPrint(int min)
virtual voidShowMembers(TMemberInspector&)
virtual voidStreamer(TBuffer&)
voidStreamerNVirtual(TBuffer& ClassDef_StreamerNVirtual_b)
RooRealVar*var(const string& varName) const
protected:
void_GetNodes(RooStats::ModelConfig* mc)
void_GetNodes(RooAbsPdf* model, const RooArgSet* observables)
map<std::string,RooAbsReal*>GetSampleFunctionMap(const string& channel)
TH1*MakeHistFromRooFunction(RooAbsReal* func, RooArgList vars, string name = "Hist")
voidPrintMultiDimHist(TH1* hist, int bin_print_width)
voidSetPrintWidths(const string& channel)
private:
RooArgSet_GetAllProducts(RooProduct* node)
RooAbsArg*findChild(const string& name, RooAbsReal* parent) const

Data Members

private:
int_bin_print_width
int_label_print_width
int_maxBinToPrint
int_minBinToPrint
vector<std::string>fChannelNameVec
map<std::string,RooArgSet*>fChannelObservMap
map<std::string,RooAbsPdf*>fChannelPdfMap
map<std::string,std::map<std::string,RooAbsReal*> >fChannelSampleFunctionMap
map<std::string,RooAbsPdf*>fChannelSumNodeMap
RooAbsPdf*fModel
RooArgSet*fObservables

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

HistFactoryNavigation(ModelConfig* mc)
HistFactoryNavigation(const string& File, const string& WorkspaceName = "combined", const string& ModelConfigName = "ModelConfig")
HistFactoryNavigation(RooAbsPdf* model, RooArgSet* observables)
void PrintMultiDimHist(TH1* hist, int bin_print_width)
RooAbsPdf* GetChannelPdf(const string& channel)
void PrintState(const string& channel)
void PrintState()
 Loop over channels and print their states, one after another
void SetPrintWidths(const string& channel)
void PrintDataSet(RooDataSet* data, const string& channel = "")
 Print the contents of a 'HistFactory' RooDataset
 These are stored in a somewhat odd way that makes
 them difficult to inspect for humans.
 They have the following layout:

 ChannelA      ChannelB     ChannelCat   Weight

 bin_1_center   0           ChannelA     bin_1_height
 bin_2_center   0           ChannelA     bin_2_height
      0        bin_1_center ChannelB     bin_1_height
      0        bin_2_center ChannelB     bin_2_height
                        ...etc...

 int label_print_width = 20;
 int bin_print_width = 12;
 Get the Data Histogram for this channel
void PrintModelAndData(RooDataSet* data)
 Loop over all channels and print model
 (including all samples) and compare
 it to the supplied dataset
void PrintParameters(bool IncludeConstantParams = false)
void PrintChannelParameters(const string& channel, bool IncludeConstantParams = false)
 Get the list of parameters
void PrintSampleParameters(const string& channel, const string& sample, bool IncludeConstantParams = false)
 Get the list of parameters
double GetBinValue(int bin, const string& channel)
 Get the total bin height for the ith bin (ROOT indexing convention)
 in channel 'channel'
 (Could be optimized, it uses an intermediate histogram for now...)
double GetBinValue(int bin, const string& channel, const string& sample)
 Get the total bin height for the ith bin (ROOT indexing convention)
 in channel 'channel'
 (This will be slow if you plan on looping over it.
  Could be optimized, it uses an intermediate histogram for now...)
std::map< std::string, RooAbsReal*> GetSampleFunctionMap(const string& channel)
 Get a map of strings to function pointers,
 which each function cooresponds to a sample
RooAbsReal* SampleFunction(const string& channel, const string& sample)
 Return the function object pointer cooresponding
 to a particular sample in a particular channel
RooArgSet* GetObservableSet(const string& channel)
 Get the observables for a particular channel
TH1* GetSampleHist(const string& channel, const string& sample, const string& name = "")
 Get a histogram of the expected values for
 a particular sample in a particular channel
 Give a name, or a default one will be used
TH1* GetChannelHist(const string& channel, const string& name = "")
 Get a histogram of the total expected value
 per bin for this channel
 Give a name, or a default one will be used
std::vector< std::string > GetChannelSampleList(const string& channel)
THStack* GetChannelStack(const string& channel, const string& name = "")
TH1* GetDataHist(RooDataSet* data, const string& channel, const string& name = "")
 TO DO:
 MAINTAIN THE ACTUAL RANGE, USING THE OBSERVABLES
 MAKE IT WORK FOR MULTI-DIMENSIONAL

 If the dataset covers multiple categories,
 Split the dataset based on the categories
void DrawChannel(const string& channel, RooDataSet* data = NULL)
RooArgSet _GetAllProducts(RooProduct* node)
void _GetNodes(RooAbsPdf* model, const RooArgSet* observables)
RooAbsArg* findChild(const string& name, RooAbsReal* parent) const
RooAbsReal* GetConstraintTerm(const string& parameter)
double GetConstraintUncertainty(const string& parameter)
void ReplaceNode(const string& ToReplace, RooAbsArg* ReplaceWith)
void PrintSampleComponents(const string& channel, const string& sample)
 Get the Sample Node
TH1* MakeHistFromRooFunction(RooAbsReal* func, RooArgList vars, string name = "Hist")
 Turn a RooAbsReal* into a TH1* based
 on a template histogram.
 The 'vars' arg list defines the x (and y and z variables)
 Loop over the bins of the Template,
 find the bin centers,
 Scan the input Var over those bin centers,
 and use the value of the function
 to make the new histogram
 Make the new histogram
 Cone and empty the template
      TH1* hist = (TH1*) histTemplate.Clone( name.c_str() );
void _GetNodes(RooStats::ModelConfig* mc)
void SetConstant(const string& regExpr = ".*", bool constant = true)
RooRealVar* var(const string& varName) const
HistFactoryNavigation(ModelConfig* mc)
 Initialze based on an already-created HistFactory Model
virtual ~HistFactoryNavigation()
{}
void SetMaxBinToPrint(int max)
{ _maxBinToPrint = max; }
int GetMaxBinToPrint() const
{ return _maxBinToPrint; }
void SetMinBinToPrint(int min)
{ _minBinToPrint = min; }
int GetMinBinToPrint() const
{ return _minBinToPrint; }
RooAbsPdf* GetModel() const
 Get the model for this channel
{ return fModel; }