A++ » ROOSTATS » 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<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& insp) const
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<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<string>fChannelNameVec
map<string,RooArgSet*>fChannelObservMap
map<string,RooAbsPdf*>fChannelPdfMap
map<string,map<string,RooAbsReal*> >fChannelSampleFunctionMap
map<string,RooAbsPdf*>fChannelSumNodeMap
RooAbsPdf*fModel
RooArgSet*fObservables

Class Charts

Inheritance Chart:
RooStats::HistFactory::HistFactoryNavigation

Function documentation

HistFactoryNavigation(ModelConfig* mc)
 Initialze based on an already-created HistFactory Model
HistFactoryNavigation(const string& File, const string& WorkspaceName = "combined", const string& ModelConfigName = "ModelConfig")
virtual ~HistFactoryNavigation()
{}
void PrintState()
 Should pretty print all channels and the current values

void PrintState(const string& channel)
 Should pretty print this and the current values
void PrintParameters(bool IncludeConstantParams = false)
 Print the current values and errors of pdf parameters
void PrintChannelParameters(const string& channel, bool IncludeConstantParams = false)
 Print parameters that effect a particular channel
void PrintSampleComponents(const string& channel, const string& sample)
 Print parameters that effect a particular sample
 Print the different components that make up a sample
 (NormFactors, Statistical Uncertainties, Interpolation, etc)
void PrintDataSet(RooDataSet* data, const string& channel = "")
 Print a "HistFactory style" RooDataSet in a readable way
void PrintModelAndData(RooDataSet* data)
 Print the model and the data, comparing channel by channel
double GetBinValue(int bin, const string& channel)
 The value of the ith bin for the total in that channel
double GetBinValue(int bin, const string& channel, const string& sample)
 The value of the ith bin for that sample and channel

TH1* GetSampleHist(const string& channel, const string& sample, const string& name = "")
 The (current) histogram for that sample
 This includes all parameters and interpolation
TH1* GetDataHist(RooDataSet* data, const string& channel, const string& name = "")
 Get the total channel histogram for this channel
 Get a histogram from the dataset for this channel
THStack* GetChannelStack(const string& channel, const string& name = "")
 Get a stack of all samples in a channel
void DrawChannel(const string& channel, RooDataSet* data = __null)
 Draw a stack of the channel, and include data if the pointer is supplied
RooAbsReal* SampleFunction(const string& channel, const string& sample)
 Get the RooAbsReal function for a given sample in a given channel
RooArgSet* GetObservableSet(const string& channel)
 Get the set of observables for a given channel
RooAbsReal* GetConstraintTerm(const string& parameter)
 Get the constraint term for a given systematic (alpha or gamma)
double GetConstraintUncertainty(const string& parameter)
 Get the uncertainty based on the constraint term for a given systematic
void ReplaceNode(const string& ToReplace, RooAbsArg* ReplaceWith)
 Find a node in the pdf and replace it with a new node
 These nodes can be functions, pdf's, RooRealVar's, etc
 Will do minimial checking to make sure the replacement makes sense
void SetConstant(const string& regExpr = ".*", bool constant = true)
 Set any RooRealVar's const (or not const) if they match
 the supplied regular expression
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; }
RooAbsPdf* GetChannelPdf(const string& channel)
std::vector< std::string > GetChannelSampleList(const string& channel)
RooRealVar* var(const string& varName) const
 Return the RooRealVar by the same name used in the model
 If not found, return NULL
void SetPrintWidths(const string& channel)
      // Add a channel to the pdf
      // Combine the data if it is provided
      void AddChannel(const std::string& channel, RooAbsPdf* pdf, RooDataSet* data=NULL);


      // Add a constraint term to the pdf
      // This method requires that the pdf NOT be simultaneous
      void AddConstraintTerm(RooAbsArg* constraintTerm);

      // Add a constraint term to the pdf of a particular channel
      // This method requires that the pdf be simultaneous
      // OR that the channel string match the channel that the pdf represents
      void AddConstraintTerm(RooAbsArg* constraintTerm, const std::string& channel);

 Set the title and bin widths
void _GetNodes(RooStats::ModelConfig* mc)
 Fetch the node information for the pdf in question, and
 save it in the varous collections in this class
void _GetNodes(RooAbsPdf* model, const RooArgSet* observables)
void PrintMultiDimHist(TH1* hist, int bin_print_width)
 Print a histogram's contents to the screen
 void PrettyPrintHistogram(TH1* hist);
TH1* MakeHistFromRooFunction(RooAbsReal* func, RooArgList vars, string name = "Hist")
 Make a histogram from a funciton
 Edit so it can take a RooArgSet of parameters
std::map< std::string, RooAbsReal*> GetSampleFunctionMap(const string& channel)
 Get a map of sample names to their functions for a particular channel
RooAbsArg* findChild(const string& name, RooAbsReal* parent) const
 Internal method implementation of finding a daughter node
 from a parent node (looping over all generations)
RooArgSet _GetAllProducts(RooProduct* node)
 Recursively get all products of products