Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ToyMCImportanceSampler.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Sven Kreiss January 2012
3// Author: Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke
4/*************************************************************************
5 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
6 * All rights reserved. *
7 * *
8 * For the licensing terms see $ROOTSYS/LICENSE. *
9 * For the list of contributors see $ROOTSYS/README/CREDITS. *
10 *************************************************************************/
11
12/** \class RooStats::ToyMCImportanceSampler
13 \ingroup Roostats
14
15ToyMCImportanceSampler is an extension of the ToyMCSampler for Importance Sampling.
16
17Implementation based on a work by Cranmer, Kreiss, Read (in Preparation)
18*/
19
21
22#include "RooMsgService.h"
23
24#include "RooCategory.h"
25#include "TMath.h"
26
27using namespace RooFit;
28using std::endl, std::vector;
29
30
31
32namespace RooStats {
33
34////////////////////////////////////////////////////////////////////////////////
35
37 for( unsigned int i=0; i < fImportanceSnapshots.size(); i++ ) if(fImportanceSnapshots[i]) delete fImportanceSnapshots[i];
38 for( unsigned int i=0; i < fNullSnapshots.size(); i++ ) if(fNullSnapshots[i]) delete fNullSnapshots[i];
39}
40
41////////////////////////////////////////////////////////////////////////////////
42
45
46 for( unsigned int i=0; i < fImpNLLs.size(); i++ ) { fImpNLLs[i].reset(); }
47 for( unsigned int i=0; i < fNullNLLs.size(); i++ ) { fNullNLLs[i].reset(); }
48}
49
50////////////////////////////////////////////////////////////////////////////////
51
53 if( fNToys == 0 ) return nullptr;
54
55 // remember original #toys, but overwrite it temporarily with the #toys per distribution
57
58 // to keep track of which dataset entry comes from which density, define a roocategory as a label
59 RooCategory densityLabel( "densityLabel", "densityLabel" );
60 densityLabel.defineType( "null", -1 );
61 for( unsigned int i=0; i < fImportanceDensities.size(); i++ )
62 densityLabel.defineType( TString::Format( "impDens_%d", i ), i );
63
64
65 RooDataSet* fullResult = nullptr;
66
67 // generate null (negative i) and imp densities (0 and positive i)
68 for( int i = -1; i < (int)fImportanceDensities.size(); i++ ) {
69 if( i < 0 ) {
70 // generate null toys
71 oocoutP(nullptr,Generation) << std::endl << std::endl << " GENERATING FROM nullptr DENSITY " << std::endl << std::endl;
72 SetDensityToGenerateFromByIndex( 0, true ); // true = generate from null
73 }else{
74 oocoutP(nullptr,Generation) << std::endl << std::endl << " GENERATING IMP DENS/SNAP "<<i+1<<" OUT OF "<<fImportanceDensities.size()<< std::endl<< std::endl;
75 SetDensityToGenerateFromByIndex( i, false ); // false = generate not from null
76 }
77
78 RooRealVar reweight( "reweight", "reweight", 1.0 );
79 // apply strategy for how to distribute the #toys between the distributions
81 // assuming alltoys = one null + N imp densities. And round up.
82 fNToys = TMath::CeilNint( double(allToys)/(fImportanceDensities.size()+1) );
84 // for N densities, split the toys into (2^(N+1))-1 parts, and assign 2^0 parts to the first
85 // density (which is the null), 2^1 to the second (first imp dens), etc, up to 2^N
86 fNToys = TMath::CeilNint( double(allToys) * pow( double(2) , i+1 ) / (pow( double(2), int(fImportanceDensities.size()+1) )-1) );
87
88 int largestNToys = TMath::CeilNint( allToys * pow( double(2), int(fImportanceDensities.size()) ) / (pow( double(2), int(fImportanceDensities.size()+1) )-1) );
89 reweight.setVal( ((double)largestNToys) / fNToys );
90 }
91
92 ooccoutI(nullptr,InputArguments) << "Generating " << fNToys << " toys for this density." << std::endl;
93 ooccoutI(nullptr,InputArguments) << "Reweight is " << reweight.getVal() << std::endl;
94
95
97
98 if (result->get()->size() > fTestStatistics.size()) {
99 // add label
100 densityLabel.setIndex( i );
101 result->addColumn( densityLabel );
102 result->addColumn( reweight );
103 }
104
105 if( !fullResult ) {
106 fullResult = new RooDataSet( result->GetName(), result->GetTitle(), *result->get(), RooFit::WeightVar());
107 }
108
109 for( int j=0; j < result->numEntries(); j++ ) {
110// std::cout << "entry: " << j << std::endl;
111// result->get(j)->Print();
112// std::cout << "weight: " << result->weight() << std::endl;
113// std::cout << "reweight: " << reweight.getVal() << std::endl;
114 const RooArgSet* row = result->get(j);
115 fullResult->add( *row, result->weight()*reweight.getVal() );
116 }
117 delete result;
118 }
119
120 // restore #toys
121 fNToys = allToys;
122
123 return fullResult;
124}
125
126////////////////////////////////////////////////////////////////////////////////
127
130 double& weight
131) const {
132 if( fNullDensities.size() > 1 ) {
133 ooccoutI(nullptr,InputArguments) << "Null Densities:" << std::endl;
134 for( unsigned int i=0; i < fNullDensities.size(); i++) {
135 ooccoutI(nullptr,InputArguments) << " null density["<<i<<"]: " << fNullDensities[i] << " \t null snapshot["<<i<<"]: " << fNullSnapshots[i] << std::endl;
136 }
137 ooccoutE(nullptr,InputArguments) << "Cannot use multiple null densities and only ask for one weight." << std::endl;
138 return nullptr;
139 }
140
141 if( fNullDensities.empty() && fPdf ) {
142 ooccoutI(nullptr,InputArguments) << "No explicit null densities specified. Going to add one based on the given paramPoint and the global fPdf. ... but cannot do that inside const function." << std::endl;
143 //AddNullDensity( fPdf, &paramPoint );
144 }
145
146 // do not do anything if the given parameter point if fNullSnapshots[0]
147 // ... which is the most common case
148 if( fNullSnapshots[0] != &paramPoint ) {
149 ooccoutD(nullptr,InputArguments) << "Using given parameter point. Replaces snapshot for the only null currently defined." << std::endl;
150 if(fNullSnapshots[0]) delete fNullSnapshots[0];
151 fNullSnapshots.clear();
152 fNullSnapshots.push_back( (RooArgSet*)paramPoint.snapshot() );
153 }
154
155 vector<double> weights;
156 weights.push_back( weight );
157
159 for( unsigned int i=0; i < fImportanceDensities.size(); i++ ) impNLLs.push_back( 0.0 );
161 for( unsigned int i=0; i < fNullDensities.size(); i++ ) nullNLLs.push_back( 0.0 );
162
164 weight = weights[0];
165 return d;
166}
167
168////////////////////////////////////////////////////////////////////////////////
169
172 double& weight,
174 double& nullNLL
175) const {
176 if( fNullDensities.size() > 1 ) {
177 ooccoutI(nullptr,InputArguments) << "Null Densities:" << std::endl;
178 for( unsigned int i=0; i < fNullDensities.size(); i++) {
179 ooccoutI(nullptr,InputArguments) << " null density["<<i<<"]: " << fNullDensities[i] << " \t null snapshot["<<i<<"]: " << fNullSnapshots[i] << std::endl;
180 }
181 ooccoutE(nullptr,InputArguments) << "Cannot use multiple null densities and only ask for one weight and NLL." << std::endl;
182 return nullptr;
183 }
184
185 if( fNullDensities.empty() && fPdf ) {
186 ooccoutI(nullptr,InputArguments) << "No explicit null densities specified. Going to add one based on the given paramPoint and the global fPdf. ... but cannot do that inside const function." << std::endl;
187 //AddNullDensity( fPdf, &paramPoint );
188 }
189
190 ooccoutI(nullptr,InputArguments) << "Using given parameter point. Overwrites snapshot for the only null currently defined." << std::endl;
191 if(fNullSnapshots[0]) delete fNullSnapshots[0];
192 fNullSnapshots.clear();
193 fNullSnapshots.push_back( (const RooArgSet*)paramPoint.snapshot() );
194
195 vector<double> weights;
196 weights.push_back( weight );
197
199 nullNLLs.push_back( nullNLL );
200
202 weight = weights[0];
203 nullNLL = nullNLLs[0];
204 return d;
205}
206
207////////////////////////////////////////////////////////////////////////////////
208
210 vector<double>& weights
211) const {
212 if( fNullDensities.size() != weights.size() ) {
213 ooccoutI(nullptr,InputArguments) << "weights.size() != nullDesnities.size(). You need to provide a vector with the correct size." << std::endl;
214 //AddNullDensity( fPdf, &paramPoint );
215 }
216
218 for( unsigned int i=0; i < fImportanceDensities.size(); i++ ) impNLLs.push_back( 0.0 );
220 for( unsigned int i=0; i < fNullDensities.size(); i++ ) nullNLLs.push_back( 0.0 );
221
223 return d;
224}
225
226////////////////////////////////////////////////////////////////////////////////
227/// This method generates a toy data set for importance sampling for the given parameter point taking
228/// global observables into account.
229/// The values of the generated global observables remain in the pdf's variables.
230/// They have to have those values for the subsequent evaluation of the
231/// test statistics.
232
234 vector<double>& weights,
237) const {
238
239
240 ooccoutD(nullptr,InputArguments) << std::endl;
241 ooccoutD(nullptr,InputArguments) << "GenerateToyDataImportanceSampling" << std::endl;
242
243 if(!fObservables) {
244 ooccoutE(nullptr,InputArguments) << "Observables not set." << std::endl;
245 return nullptr;
246 }
247
248 if( fNullDensities.empty() ) {
249 oocoutE(nullptr,InputArguments) << "ToyMCImportanceSampler: Need to specify the null density explicitly." << std::endl;
250 return nullptr;
251 }
252
253 // catch the case when NLLs are not created (e.g. when ToyMCSampler was streamed for Proof)
254 if( fNullNLLs.empty() && !fNullDensities.empty() ) {
255 for( unsigned int i = 0; i < fNullDensities.size(); i++ ) fNullNLLs.push_back( nullptr );
256 }
257 if( fImpNLLs.empty() && !fImportanceDensities.empty() ) {
258 for( unsigned int i = 0; i < fImportanceDensities.size(); i++ ) fImpNLLs.push_back( nullptr );
259 }
260
261 if( fNullDensities.size() != fNullNLLs.size() ) {
262 oocoutE(nullptr,InputArguments) << "ToyMCImportanceSampler: Something wrong. NullNLLs must be of same size as null densities." << std::endl;
263 return nullptr;
264 }
265
268 ) {
269 oocoutE(nullptr,InputArguments) << "ToyMCImportanceSampler: no importance density given or index out of range." << std::endl;
270 return nullptr;
271 }
272
273
274 // paramPoint used to be given as parameter
275 // situation is clear when there is only one null.
276 // WHAT TO DO FOR MANY nullptr DENSITIES?
278 //cout << "paramPoint: " << std::endl;
279 //paramPoint.Print("v");
280
281
282 // assign input paramPoint
283 std::unique_ptr<RooArgSet> allVars{fPdf->getVariables()};
284 allVars->assign(paramPoint);
285
286
287 // create nuisance parameter points
290
291 // generate global observables
292 RooArgSet observables(*fObservables);
293 if(fGlobalObservables && !fGlobalObservables->empty()) {
294 observables.remove(*fGlobalObservables);
295 // WHAT TODO FOR MANY nullptr DENSITIES?
297 }
298
299 // save values to restore later.
300 // but this must remain after(!) generating global observables
301 if( !fGenerateFromNull ) {
302 std::unique_ptr<RooArgSet> allVarsImpDens{fImportanceDensities[fIndexGenDensity]->getVariables()};
303 allVars->add(*allVarsImpDens);
304 }
305 const RooArgSet* saveVars = (const RooArgSet*)allVars->snapshot();
306
307 double globalWeight = 1.0;
308 if(fNuisanceParametersSampler) { // use nuisance parameters?
309 // Construct a set of nuisance parameters that has the parameters
310 // in the input paramPoint removed. Therefore, no parameter in
311 // paramPoint is randomized.
312 // Therefore when a parameter is given (should be held fixed),
313 // but is also in the list of nuisance parameters, the parameter
314 // will be held fixed. This is useful for debugging to hold single
315 // parameters fixed although under "normal" circumstances it is
316 // randomized.
318 allVarsMinusParamPoint.remove(paramPoint, false, true); // match by name
319
320 // get nuisance parameter point and weight
322 }
323 // populate input weights vector with this globalWeight
324 for( unsigned int i=0; i < weights.size(); i++ ) weights[i] = globalWeight;
325
326 RooAbsData* data = nullptr;
327 if( fGenerateFromNull ) {
328 //cout << "gen from null" << std::endl;
329 allVars->assign(*fNullSnapshots[fIndexGenDensity]);
330 data = Generate(*fNullDensities[fIndexGenDensity], observables).release();
331 }else{
332 // need to be careful here not to overwrite the current state of the
333 // nuisance parameters, ie they must not be part of the snapshot
334 //cout << "gen from imp" << std::endl;
336 allVars->assign(*fImportanceSnapshots[fIndexGenDensity]);
337 }
338 data = Generate(*fImportanceDensities[fIndexGenDensity], observables).release();
339 }
340 //cout << "data generated: " << data << std::endl;
341
342 if (!data) {
343 oocoutE(nullptr,InputArguments) << "ToyMCImportanceSampler: error generating data" << std::endl;
344 return nullptr;
345 }
346
347
348
349 // Importance Sampling: adjust weight
350 // Sources: Alex Read, presentation by Michael Woodroofe
351
352 ooccoutD(nullptr,InputArguments) << "About to create/calculate all nullNLLs." << std::endl;
353 for( unsigned int i=0; i < fNullDensities.size(); i++ ) {
354 //oocoutI(nullptr,InputArguments) << "Setting variables to nullSnapshot["<<i<<"]"<< std::endl;
355 //fNullSnapshots[i]->Print("v");
356
357 allVars->assign(*fNullSnapshots[i]);
358 if( !fNullNLLs[i] ) {
359 std::unique_ptr<RooArgSet> allParams{fNullDensities[i]->getParameters(*data)};
360 fNullNLLs[i] = std::unique_ptr<RooAbsReal>{fNullDensities[i]->createNLL(*data, RooFit::CloneData(false), RooFit::Constrain(*allParams),
362 }else{
363 fNullNLLs[i]->setData( *data, false );
364 }
365 nullNLLVals[i] = fNullNLLs[i]->getVal();
366 // FOR DEBuGGING!!!!!!!!!!!!!!!!!
367 if( !fReuseNLL ) { fNullNLLs[i].reset(); }
368 }
369
370
371 // for each null: find minNLLVal of null and all imp densities
372 ooccoutD(nullptr,InputArguments) << "About to find the minimum NLLs." << std::endl;
374 for( unsigned int i=0; i < nullNLLVals.size(); i++ ) minNLLVals.push_back( nullNLLVals[i] );
375
376 for( unsigned int i=0; i < fImportanceDensities.size(); i++ ) {
377 //oocoutI(nullptr,InputArguments) << "Setting variables to impSnapshot["<<i<<"]"<< std::endl;
378 //fImportanceSnapshots[i]->Print("v");
379
380 if( fImportanceSnapshots[i] ) {
381 allVars->assign(*fImportanceSnapshots[i]);
382 }
383 if( !fImpNLLs[i] ) {
384 std::unique_ptr<RooArgSet> allParams{fImportanceDensities[i]->getParameters(*data)};
385 fImpNLLs[i] = std::unique_ptr<RooAbsReal>{fImportanceDensities[i]->createNLL(*data, RooFit::CloneData(false), RooFit::Constrain(*allParams),
387 }else{
388 fImpNLLs[i]->setData( *data, false );
389 }
390 impNLLVals[i] = fImpNLLs[i]->getVal();
391 // FOR DEBuGGING!!!!!!!!!!!!!!!!!
392 if( !fReuseNLL ) { fImpNLLs[i].reset(); }
393
394 for( unsigned int j=0; j < nullNLLVals.size(); j++ ) {
395 if( impNLLVals[i] < minNLLVals[j] ) minNLLVals[j] = impNLLVals[i];
396 ooccoutD(nullptr,InputArguments) << "minNLLVals["<<j<<"]: " << minNLLVals[j] << " nullNLLVals["<<j<<"]: " << nullNLLVals[j] << " impNLLVals["<<i<<"]: " << impNLLVals[i] << std::endl;
397 }
398 }
399
400 // veto toys: this is a sort of "overlap removal" of the various distributions
401 // if not vetoed: apply weight
402 ooccoutD(nullptr,InputArguments) << "About to apply vetos and calculate weights." << std::endl;
403 for( unsigned int j=0; j < nullNLLVals.size(); j++ ) {
404 if ( fApplyVeto && fGenerateFromNull && minNLLVals[j] != nullNLLVals[j] ) weights[j] = 0.0;
405 else if( fApplyVeto && !fGenerateFromNull && minNLLVals[j] != impNLLVals[fIndexGenDensity] ) weights[j] = 0.0;
406 else if( !fGenerateFromNull ) {
407 // apply (for fImportanceGenNorm, the weight is one, so nothing needs to be done here)
408
409 // L(pdf) / L(imp) = exp( NLL(imp) - NLL(pdf) )
410 weights[j] *= exp(minNLLVals[j] - nullNLLVals[j]);
411 }
412
413 ooccoutD(nullptr,InputArguments) << "weights["<<j<<"]: " << weights[j] << std::endl;
414 }
415
416
417
418 allVars->assign(*saveVars);
419 delete saveVars;
420
421 return data;
422}
423
424////////////////////////////////////////////////////////////////////////////////
425/// poi has to be fitted beforehand. This function expects this to be the muhat value.
426
428 // these might not necessarily be the same thing.
429 double impMaxMu = poi.getVal();
430
431 // this includes the null
432 int n = 1;
433
434 // check whether error is trustworthy
435 if( poi.getError() > 0.01 && poi.getError() < 5.0 ) {
436 n = TMath::CeilNint( poi.getVal() / (2.*nStdDevOverlap*poi.getError()) ); // round up
437 oocoutI(nullptr,InputArguments) << "Using fitFavoredMu and error to set the number of imp points" << std::endl;
438 oocoutI(nullptr,InputArguments) << "muhat: " << poi.getVal() << " optimize for distance: " << 2.*nStdDevOverlap*poi.getError() << std::endl;
439 oocoutI(nullptr,InputArguments) << "n = " << n << std::endl;
440 oocoutI(nullptr,InputArguments) << "This results in a distance of: " << impMaxMu / n << std::endl;
441 }
442
443 // exclude the null, just return the number of importance snapshots
445}
446
447////////////////////////////////////////////////////////////////////////////////
448/// n is the number of importance densities
449
451
452 // these might not necessarily be the same thing.
453 double impMaxMu = poi.getVal();
454
455 // create imp snapshots
456 if( impMaxMu > poiValueForBackground && n > 0 ) {
457 for( int i=1; i <= n; i++ ) {
459 oocoutI(nullptr,InputArguments) << std::endl << "create point with poi: " << std::endl;
460 poi.Print();
461
462 // impSnaps without first snapshot because that is null hypothesis
463
465 }
466 }
467
468 return n;
469}
470
471} // end namespace RooStats
#define d(i)
Definition RSha256.hxx:102
#define oocoutE(o, a)
#define oocoutI(o, a)
#define ooccoutI(o, a)
#define ooccoutD(o, a)
#define oocoutP(o, a)
#define ooccoutE(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
void Print(Option_t *options=nullptr) const override
Print the object to the defaultPrintStream().
Definition RooAbsArg.h:263
RooFit::OwningPtr< RooArgSet > getVariables(bool stripDisconnected=true) const
Return RooArgSet with all variables (tree leaf nodes of expression tree)
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Object to represent discrete states.
Definition RooCategory.h:28
Container class to hold unbinned data.
Definition RooDataSet.h:34
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
double getError() const
Definition RooRealVar.h:58
Helper class for ToyMCSampler.
void NextPoint(RooArgSet &nuisPoint, double &weight)
Assigns new nuisance parameter point to members of nuisPoint.
void AddImportanceDensity(RooAbsPdf *p, const RooArgSet *s)
For importance sampling with multiple densities/snapshots: This is used to check the current Likeliho...
std::vector< const RooArgSet * > fImportanceSnapshots
RooDataSet * GetSamplingDistributionsSingleWorker(RooArgSet &paramPoint) override
overwrite GetSamplingDistributionsSingleWorker(paramPoint) with a version that loops over nulls and i...
RooAbsData * GenerateToyData(RooArgSet &paramPoint, double &weight) const override
std::vector< const RooArgSet * > fNullSnapshots
void SetDensityToGenerateFromByIndex(unsigned int i, bool fromNull=false)
specifies the pdf to sample from
std::vector< RooAbsPdf * > fImportanceDensities
std::vector< std::unique_ptr< RooAbsReal > > fImpNLLs
!
std::vector< std::unique_ptr< RooAbsReal > > fNullNLLs
!
int CreateImpDensitiesForOnePOIAdaptively(RooAbsPdf &pdf, const RooArgSet &allPOI, RooRealVar &poi, double nStdDevOverlap=0.5, double poiValueForBackground=0.0)
poi has to be fitted beforehand. This function expects this to be the muhat value.
void ClearCache() override
helper method for clearing the cache
RooArgSet fConditionalObs
set of conditional observables
std::vector< RooAbsPdf * > fNullDensities
support multiple null densities
int CreateNImpDensitiesForOnePOI(RooAbsPdf &pdf, const RooArgSet &allPOI, RooRealVar &poi, int n, double poiValueForBackground=0.0)
n is the number of importance densities
const RooArgSet * fGlobalObservables
virtual void GenerateGlobalObservables(RooAbsPdf &pdf) const
generate global observables
virtual RooDataSet * GetSamplingDistributionsSingleWorker(RooArgSet &paramPoint)
This is the main function for serial runs.
NuisanceParametersSampler * fNuisanceParametersSampler
!
const RooArgSet * fObservables
RooAbsPdf * fPdf
densities, snapshots, and test statistics to reweight to
bool fExpectedNuisancePar
whether to use expectation values for nuisance parameters (ie Asimov data set)
const RooArgSet * fNuisancePars
std::vector< TestStatistic * > fTestStatistics
std::unique_ptr< RooAbsData > Generate(RooAbsPdf &pdf, RooArgSet &observables, const RooDataSet *protoData=nullptr, int forceEvents=0) const
helper for GenerateToyData
Int_t fNToys
number of toys to generate
virtual void ClearCache()
helper method for clearing the cache
RooAbsPdf * fPriorNuisance
prior pdf for nuisance parameters
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
RooCmdArg WeightVar(const char *name="weight", bool reinterpretAsWeight=false)
RooCmdArg Constrain(const RooArgSet &params)
RooCmdArg CloneData(bool flag)
RooCmdArg ConditionalObservables(Args_t &&... argsOrArgSet)
Create a RooCmdArg to declare conditional observables.
const Int_t n
Definition legend1.C:16
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
@ InputArguments
Namespace for the RooStats classes.
Definition CodegenImpl.h:58
Int_t CeilNint(Double_t x)
Returns the nearest integer of TMath::Ceil(x).
Definition TMath.h:678