Logo ROOT   master
Reference Guide
RooUniform.h
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * File: $Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 #ifndef ROO_UNIFORM
17 #define ROO_UNIFORM
18 
19 #include "RooAbsPdf.h"
20 #include "RooListProxy.h"
21 
22 class RooRealVar;
23 
24 class RooUniform : public RooAbsPdf {
25 public:
26  RooUniform() {} ;
27  RooUniform(const char *name, const char *title, const RooArgSet& _x);
28  RooUniform(const RooUniform& other, const char* name=0) ;
29  virtual TObject* clone(const char* newname) const { return new RooUniform(*this,newname); }
30  inline virtual ~RooUniform() { }
31 
32  Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
33  Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
34 
35  Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
36  void generateEvent(Int_t code);
37 
38 protected:
39 
41 
42  Double_t evaluate() const ;
43  inline RooSpan<double> evaluateBatch(std::size_t, std::size_t ) const {
44  return {};
45  }
46 
47 private:
48 
49  ClassDef(RooUniform,1) // Flat PDF in N dimensions
50 };
51 
52 #endif
void generateEvent(Int_t code)
Implement internal generator.
Definition: RooUniform.cxx:127
A simple container to hold a batch of data values.
Definition: RooSpan.h:32
RooSpan< double > evaluateBatch(std::size_t, std::size_t) const
Evaluate function for a batch of input data points.
Definition: RooUniform.h:43
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Advertise analytical integral.
Definition: RooUniform.cxx:65
virtual ~RooUniform()
Definition: RooUniform.h:30
#define ClassDef(name, id)
Definition: Rtypes.h:326
RooRealVar represents a fundamental (non-derived) real-valued object.
Definition: RooRealVar.h:36
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implement analytical integral.
Definition: RooUniform.cxx:88
RooListProxy x
Definition: RooUniform.h:40
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooUniform.cxx:57
Flat p.d.f.
Definition: RooUniform.h:24
RooListProxy is the concrete proxy for RooArgList objects.
Definition: RooListProxy.h:25
Mother of all ROOT objects.
Definition: TObject.h:37
RooAbsPdf, the base class of all PDFs
Definition: RooAbsPdf.h:40
const Bool_t kTRUE
Definition: RtypesCore.h:87
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Advertise internal generator.
Definition: RooUniform.cxx:103
char name[80]
Definition: TGX11.cxx:109
virtual TObject * clone(const char *newname) const
Definition: RooUniform.h:29