+
class RooMCStudy
-
library: libRooFitCore
#include "RooMCStudy.h"
Display options:
Show inherited
Show non-public

class RooMCStudy

Function Members (Methods)

public:
RooMCStudy(const RooAbsPdf& genModel, const RooAbsPdf& fitModel, const RooArgSet& dependents, const char* genOptions = "", const char* fitOptions = "", const RooDataSet* genProtoData = 0, const RooArgSet& projDeps = RooArgSet())
RooMCStudy(const RooAbsPdf& model, const RooArgSet& observables, RooCmdArg arg1 = RooCmdArg::none, RooCmdArg arg2 = RooCmdArg::none, RooCmdArg arg3 = RooCmdArg::none, RooCmdArg arg4 = RooCmdArg::none, RooCmdArg arg5 = RooCmdArg::none, RooCmdArg arg6 = RooCmdArg::none, RooCmdArg arg7 = RooCmdArg::none, RooCmdArg arg8 = RooCmdArg::none)
virtual~RooMCStudy()
Bool_taddFitResult(const RooFitResult& fr)
voidaddModule(RooAbsMCStudyModule& module)
static TClass*Class()
Bool_tfit(Int_t nSamples, const char* asciiFilePat)
Bool_tfit(Int_t nSamples, TList& dataSetList)
const RooArgSet*fitParams(Int_t sampleNum) const
const RooDataSet&fitParDataSet()
const RooFitResult*fitResult(Int_t sampleNum) const
const RooDataSet*genData(Int_t sampleNum) const
Bool_tgenerate(Int_t nSamples, Int_t nEvtPerSample = 0, Bool_t keepGenData = kFALSE, const char* asciiFilePat = 0)
Bool_tgenerateAndFit(Int_t nSamples, Int_t nEvtPerSample = 0, Bool_t keepGenData = kFALSE, const char* asciiFilePat = 0)
virtual TClass*IsA() const
RooPlot*plotError(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins = 100)
RooPlot*plotError(const RooRealVar& param, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPlot*plotNLL(Double_t lo, Double_t hi, Int_t nBins = 100)
RooPlot*plotNLL(const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPlot*plotParam(const RooRealVar& param, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPlot*plotParam(const char* paramName, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPlot*plotParamOn(RooPlot* frame, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
RooPlot*plotPull(const RooRealVar& param, Double_t lo = -3.0, Double_t hi = 3.0, Int_t nbins = 25, Bool_t fitGauss = kFALSE)
RooPlot*plotPull(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
virtual voidShowMembers(TMemberInspector& insp, char* parent)
virtual voidStreamer(TBuffer& b)
voidStreamerNVirtual(TBuffer& b)
protected:
voidcalcPulls()
RooFitResult*doFit(RooAbsData* genSample)
Bool_tfitSample(RooAbsData* genSample)
RooPlot*makeFrameAndPlotCmd(const RooRealVar& param, RooLinkedList& cmdList, Bool_t symRange = kFALSE) const
RooFitResult*refit(RooAbsData* genSample = 0)
voidresetFitParams()
Bool_trun(Bool_t generate, Bool_t fit, Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
private:
RooMCStudy(const RooMCStudy&)

Data Members

protected:
RooArgSet_allDependentsList of generate + prototype dependents
Bool_t_binGenDataBin data between generating and fitting
Bool_t_canAddFitResultsAllow adding of external fit results?
RooArgSet_dependentsList of dependents
Bool_t_extendedGenAdd poisson term to number of events to generate?
RooArgSet*_fitInitParamsList of initial values of fit parameters
RooAbsPdf*_fitModelFit model
RooLinkedList_fitOptListFit option command list
TString_fitOptionsFit options string
RooDataSet*_fitParDataData set of fit parameters of each sample
RooArgSet*_fitParamsList of actual fit parameters
TList_fitResListList of RooFitResult fit output objects
RooAbsGenContext*_genContextGenerator context
TList_genDataListList of generated data sample
RooArgSet*_genInitParamsList of originalgenerator parameters
RooAbsPdf*_genModelGenerator model
RooArgSet*_genParamsList of actual generator parameters
const RooDataSet*_genProtoDataGenerator prototype data set
RooDataSet*_genSampleCurrently generated sample
list<RooAbsMCStudyModule*>_modListList of additional study modules ;
Double_t_nExpGenNumber of expected events to generate in extended mode
RooRealVar*_nllVar
RooArgSet_projDepsList of projected dependents in fit
Bool_t_randProtoRandomize order of prototype data access
Bool_t_verboseGenVerbose generation?

Class Charts

Inheritance Inherited Members Includes Libraries
Class Charts

Function documentation

RooMCStudy(const RooAbsPdf& model, const RooArgSet& observables, RooCmdArg arg1 = RooCmdArg::none, RooCmdArg arg2 = RooCmdArg::none, RooCmdArg arg3 = RooCmdArg::none, RooCmdArg arg4 = RooCmdArg::none, RooCmdArg arg5 = RooCmdArg::none, RooCmdArg arg6 = RooCmdArg::none, RooCmdArg arg7 = RooCmdArg::none, RooCmdArg arg8 = RooCmdArg::none)
 Construct Monte Carlo Study Manager. This class automates generating data from a given PDF,
 fitting the PDF to that data and accumulating the fit statistics.

 The constructor accepts the following arguments

 model       -- The PDF to be studied
 observables -- The variables of the PDF to be considered the observables

 FitModel(const RooAbsPdf&)        -- The PDF for fitting, if it is different from the PDF for generating
 ConditionalObservables
           (const RooArgSet& set)  -- The set of observables that the PDF should _not_ be normalized over
 Binned(Bool_t flag)               -- Bin the dataset before fitting it. Speeds up fitting of large data samples
 FitOptions(const char*)           -- Classic fit options, provided for backward compatibility
 FitOptions(....)                  -- Options to be used for fitting. All named arguments inside FitOptions()
                                                   are passed to RooAbsPdf::fitTo();
 Verbose(Bool_t flag)              -- Activate informational messages in event generation phase
 Extended(Bool_t flag)             -- Determine number of events for each sample anew from a Poisson distribution
 ProtoData(const RooDataSet&,
                 Bool_t randOrder) -- Prototype data for the event generation. If the randOrder flag is
                                      set, the order of the dataset will be re-randomized for each generation
                                      cycle to protect against systematic biases if the number of generated
                                      events does not exactly match the number of events in the prototype dataset
                                      at the cost of reduced precision
                                      with mu equal to the specified number of events
 Stuff all arguments in a list
RooMCStudy(const RooAbsPdf& genModel, const RooAbsPdf& fitModel, const RooArgSet& dependents, const char* genOptions = "", const char* fitOptions = "", const RooDataSet* genProtoData = 0, const RooArgSet& projDeps = RooArgSet())
 Constructor with a generator and fit model. Both models may point
 to the same object. The 'dependents' set of variables is generated
 in the generator phase. The optional prototype dataset is passed to
 the generator

 Available generator options
  v  - Verbose
  e  - Extended: use Poisson distribution for Nevts generated

 Available fit options
  See RooAbsPdf::fitTo()

~RooMCStudy()
 Destructor
void addModule(RooAbsMCStudyModule& module)
 Method to add study modules
Bool_t run(Bool_t generate, Bool_t fit, Int_t nSamples, Int_t nEvtPerSample, Bool_t keepGenData, const char* asciiFilePat)
 Run engine. Generate and/or fit, according to flags, 'nSamples' samples of 'nEvtPerSample' events.
 If keepGenData is set, all generated data sets will be kept in memory and can be accessed
 later via genData().

 When generating, data sets will be written out in ascii form if the pattern string is supplied
 The pattern, which is a template for sprintf, should look something like "data/toymc_%04d.dat"
 and should contain one integer field that encodes the sample serial number.

 When fitting only, data sets may optionally be read from ascii files, using the same file
 pattern.

Bool_t generateAndFit(Int_t nSamples, Int_t nEvtPerSample = 0, Bool_t keepGenData = kFALSE, const char* asciiFilePat = 0)
 Generate and fit 'nSamples' samples of 'nEvtPerSample' events.
 If keepGenData is set, all generated data sets will be kept in memory and can be accessed
 later via genData().

 Data sets will be written out is ascii form if the pattern string is supplied.
 The pattern, which is a template for sprintf, should look something like "data/toymc_%04d.dat"
 and should contain one integer field that encodes the sample serial number.

Bool_t generate(Int_t nSamples, Int_t nEvtPerSample = 0, Bool_t keepGenData = kFALSE, const char* asciiFilePat = 0)
 Generate 'nSamples' samples of 'nEvtPerSample' events.
 If keepGenData is set, all generated data sets will be kept in memory
 and can be accessed later via genData().

 Data sets will be written out in ascii form if the pattern string is supplied.
 The pattern, which is a template for sprintf, should look something like "data/toymc_%04d.dat"
 and should contain one integer field that encodes the sample serial number.

Bool_t fit(Int_t nSamples, const char* asciiFilePat)
 Fit 'nSamples' datasets, which are read from ASCII files.

 The ascii file pattern, which is a template for sprintf, should look something like "data/toymc_%04d.dat"
 and should contain one integer field that encodes the sample serial number.

Bool_t fit(Int_t nSamples, TList& dataSetList)
 Fit 'nSamples' datasets, as supplied in 'dataSetList'

void resetFitParams()
RooFitResult* doFit(RooAbsData* genSample)
 Perform actual fit according to specifications
RooFitResult* refit(RooAbsData* genSample = 0)
Bool_t fitSample(RooAbsData* genSample)
 Fit given dataset with fit model. If fit
 converges (TMinuit status code zero)
 The fit results are appended to the fit results
 dataset

 If the fit option "r" is supplied, the RooFitResult
 objects will always be saved, regardless of the
 fit status. RooFitResults objects can be retrieved
 later via fitResult().

Bool_t addFitResult(const RooFitResult& fr)
void calcPulls()
 Calculate the pulls for all fit parameters in
 the fit results data set, and add them to that dataset
const RooDataSet& fitParDataSet()
 Return the fit parameter dataset
const RooArgSet* fitParams(Int_t sampleNum) const
 Return an argset with the fit parameters for the given sample number
 NB: The fit parameters are only stored for successfull fits,
     thus the maximum sampleNum can be less that the number
     of generated samples and if so, the indeces will
     be out of synch with genData() and fitResult()
const RooFitResult* fitResult(Int_t sampleNum) const
 Return the fit result object of the fit to given sample
const RooDataSet* genData(Int_t sampleNum) const
 Return the given generated dataset
RooPlot* plotParamOn(RooPlot* frame, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
 Plot the distribution of the fitted value of the given parameter on the specified frame
 Any specified named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* plotParam(const char* paramName, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
 Plot the distribution of the fitted value of the given parameter on a newly created frame.

 This function accepts the following optional arguments
 FrameRange(double lo, double hi) -- Set range of frame to given specification
 FrameBins(int bins)              -- Set default number of bins of frame to given number
 Frame(...)                       -- Pass supplied named arguments to RooAbsRealLValue::frame() function. See frame() function
                                     for list of allowed arguments

 If no frame specifications are given, the AutoRange() feature will be used to set the range
 Any other named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* plotParam(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2, const RooCmdArg& arg3, const RooCmdArg& arg4, const RooCmdArg& arg5, const RooCmdArg& arg6, const RooCmdArg& arg7, const RooCmdArg& arg8)
 Plot the distribution of the fitted value of the given parameter on a newly created frame.

 This function accepts the following optional arguments
 FrameRange(double lo, double hi) -- Set range of frame to given specification
 FrameBins(int bins)              -- Set default number of bins of frame to given number
 Frame(...)                       -- Pass supplied named arguments to RooAbsRealLValue::frame() function. See frame() function
                                     for list of allowed arguments

 If no frame specifications are given, the AutoRange() feature will be used to set the range
 Any other named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* plotNLL(const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
 Plot the distribution of the -log(l) values on a newly created frame.

 This function accepts the following optional arguments
 FrameRange(double lo, double hi) -- Set range of frame to given specification
 FrameBins(int bins)              -- Set default number of bins of frame to given number
 Frame(...)                       -- Pass supplied named arguments to RooAbsRealLValue::frame() function. See frame() function
                                     for list of allowed arguments

 If no frame specifications are given, the AutoRange() feature will be used to set the range
 Any other named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* plotError(const RooRealVar& param, const RooCmdArg& arg1 = RooCmdArg::none, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
 Plot the distribution of the fit errors for the specified parameter on a newly created frame.

 This function accepts the following optional arguments
 FrameRange(double lo, double hi) -- Set range of frame to given specification
 FrameBins(int bins)              -- Set default number of bins of frame to given number
 Frame(...)                       -- Pass supplied named arguments to RooAbsRealLValue::frame() function. See frame() function
                                     for list of allowed arguments

 If no frame specifications are given, the AutoRange() feature will be used to set the range
 Any other named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* plotPull(const RooRealVar& param, const RooCmdArg& arg1, const RooCmdArg& arg2 = RooCmdArg::none, const RooCmdArg& arg3 = RooCmdArg::none, const RooCmdArg& arg4 = RooCmdArg::none, const RooCmdArg& arg5 = RooCmdArg::none, const RooCmdArg& arg6 = RooCmdArg::none, const RooCmdArg& arg7 = RooCmdArg::none, const RooCmdArg& arg8 = RooCmdArg::none)
 Plot the distribution of pull values for the specified parameter on a newly created frame. If asymmetric
 errors are calculated in the fit (by MINOS) those will be used in the pull calculation

 This function accepts the following optional arguments
 FrameRange(double lo, double hi) -- Set range of frame to given specification
 FrameBins(int bins)              -- Set default number of bins of frame to given number
 Frame(...)                       -- Pass supplied named arguments to RooAbsRealLValue::frame() function. See frame() function
                                     for list of allowed arguments
 FitGauss(Bool_t flag)            -- Add a gaussian fit to the frame

 If no frame specifications are given, the AutoSymRange() feature will be used to set the range
 Any other named argument is passed to the RooAbsData::plotOn() call. See that function for allowed options
RooPlot* makeFrameAndPlotCmd(const RooRealVar& param, RooLinkedList& cmdList, Bool_t symRange = kFALSE) const
RooPlot* plotNLL(Double_t lo, Double_t hi, Int_t nBins = 100)
 Create a RooPlot of the NLL distribution in the range lo-hi
 with 'nBins' bins
RooPlot* plotError(const RooRealVar& param, Double_t lo, Double_t hi, Int_t nbins = 100)
 Create a RooPlot of the distribution of the fitted errors of the given parameter.
 The range lo-hi is plotted in nbins bins
RooPlot* plotPull(const RooRealVar& param, Double_t lo = -3.0, Double_t hi = 3.0, Int_t nbins = 25, Bool_t fitGauss = kFALSE)
 Create a RooPlot of the pull distribution for the given parameter.
 The range lo-hi is plotted in nbins.
 If fitGauss is set, an unbinned max. likelihood fit of the distribution to a Gaussian model
 is performed. The fit result is overlaid on the returned RooPlot and a box with the fitted
 mean and sigma is added.
RooMCStudy(const RooAbsPdf& model, const RooArgSet& observables, RooCmdArg arg1 = RooCmdArg::none, RooCmdArg arg2 = RooCmdArg::none, RooCmdArg arg3 = RooCmdArg::none, RooCmdArg arg4 = RooCmdArg::none, RooCmdArg arg5 = RooCmdArg::none, RooCmdArg arg6 = RooCmdArg::none, RooCmdArg arg7 = RooCmdArg::none, RooCmdArg arg8 = RooCmdArg::none)

Last update: Mon Jun 25 19:44:32 2007
Copyright (c) 2000-2005, Regents of the University of California *

This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.