Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf902_numgenconfig.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Numeric algorithm tuning: configuration and customization of how MC sampling algorithms on specific pdfs are executed

#include "RooRealVar.h"
#include "RooDataSet.h"
#include "RooChebychev.h"
#include "TCanvas.h"
#include "TAxis.h"
#include "RooPlot.h"
#include "RooArgSet.h"
#include <iomanip>
using namespace RooFit;
{
// A d j u s t g l o b a l MC s a m p l i n g s t r a t e g y
// ------------------------------------------------------------------
// Example pdf for use below
RooRealVar x("x", "x", 0, 10);
RooChebychev model("model", "model", x, RooArgList(0, 0.5, -0.1));
// Change global strategy for 1D sampling problems without conditional observable
// (1st false) and without discrete observable (2nd false) from RooFoamGenerator,
// ( an interface to the TFoam MC generator with adaptive subdivisioning strategy ) to RooAcceptReject,
// a plain accept/reject sampling algorithm [ RooFit default before ROOT 5.23/04 ]
RooAbsPdf::defaultGeneratorConfig()->method1D(false, false).setLabel("RooAcceptReject");
// Generate 10Kevt using RooAcceptReject
std::unique_ptr<RooDataSet> data_ar{model.generate(x, 10000, Verbose(true))};
data_ar->Print();
// A d j u s t i n g d e f a u l t c o n f i g f o r a s p e c i f i c p d f
// -------------------------------------------------------------------------------------
// Another possibility: associate custom MC sampling configuration as default for object 'model'
// The true argument will install a clone of the default configuration as specialized configuration
// for this model if none existed so far
model.specialGeneratorConfig(true)->method1D(false, false).setLabel("RooFoamGenerator");
// A d j u s t i n g p a r a m e t e r s o f a s p e c i f i c t e c h n i q u e
// ---------------------------------------------------------------------------------------
// Adjust maximum number of steps of RooIntegrator1D in the global default configuration
RooAbsPdf::defaultGeneratorConfig()->getConfigSection("RooAcceptReject").setRealValue("nTrial1D", 2000);
// Example of how to change the parameters of a numeric integrator
// (Each config section is a RooArgSet with RooRealVars holding real-valued parameters
// and RooCategories holding parameters with a finite set of options)
model.specialGeneratorConfig()->getConfigSection("RooFoamGenerator").setRealValue("chatLevel", 1);
// Generate 10Kevt using RooFoamGenerator (FOAM verbosity increased with above chatLevel adjustment for illustration
// purposes)
std::unique_ptr<RooDataSet> data_foam{model.generate(x, 10000, Verbose())};
data_foam->Print();
}
bool setRealValue(const char *name, double newVal=0.0, bool verbose=false)
Set value of a RooAbsRealLValue stored in set with given name to newVal No error messages are printed...
static RooNumGenConfig * defaultGeneratorConfig()
Returns the default numeric MC generator configuration for all RooAbsReals.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
bool setLabel(const char *label, bool printError=true) override
Set value by specifying the name of the desired state.
Chebychev polynomial p.d.f.
RooCategory & method1D(bool cond, bool cat)
const RooArgSet & getConfigSection(const char *name) const
Retrieve configuration information specific to integrator with given name.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Double_t x[n]
Definition legend1.C:17
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
[#0] ERROR:InputArguments -- Trying to set invalid state label 'RooAcceptReject' for category method1D
--- RooGenContext ---
Using PDF RooChebychev::model[ x=x coefList=(0,0.5,-0.1) ]
Use PDF generator for ()
Use MC sampling generator RooFoamGenerator for (x)
RooDataSet::modelData[x] = 10000 entries
FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF
F F
F **************************************** F
F ****** TFoam::Initialize ****** F
F **************************************** F
F TFOAM F
F Version = 1.02M = Release date: 2005.04.10 F
F kDim = 1 = Dimension of the hyper-cubical space F
F nCells = 30 = Requested number of Cells (half of them active) F
F nSampl = 200 = No of MC events in exploration of a cell F
F nBin = 8 = No of bins in histograms, MC exploration of cell F
F EvPerBin = 25 = Maximum No effective_events/bin, MC exploration F
F OptDrive = 2 = Type of Driver =1,2 for Sigma,WtMax F
F OptRej = 1 = MC rejection on/off for OptRej=0,1 F
F MaxWtRej = 1.1 = Maximum wt in rejection for wt=1 evts F
F F
FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF
11
FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF
F F
F *** TFoam::Initialize FINISHED!!! *** F
F nCalls = 5800 = Total number of function calls F
F XPrime = 0.10992972 = Primary total integral F
F XDiver = 0.010000374 = Driver total integral F
F mcResult = 0.099929343 = Estimate of the true MC Integral F
F F
FFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFFF
--- RooGenContext ---
Using PDF RooChebychev::model[ x=x coefList=(0,0.5,-0.1) ]
Use PDF generator for ()
Use MC sampling generator RooFoamGenerator for (x)
RooDataSet::modelData[x] = 10000 entries
Date
July 2008
Author
Wouter Verkerke

Definition in file rf902_numgenconfig.C.