Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooLandau.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooLandau.h,v 1.5 2007/07/12 20:30:49 wouter Exp $
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_LANDAU
17#define ROO_LANDAU
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21
22class RooRealVar;
23
24class RooLandau : public RooAbsPdf {
25public:
26 RooLandau() {} ;
27 // Original constructor without RooAbsReal::Ref for backwards compatibility.
28 inline RooLandau(const char *name, const char *title, RooAbsReal& _x, RooAbsReal& _mean, RooAbsReal& _sigma)
29 : RooLandau{name, title, RooAbsReal::Ref{_x}, RooAbsReal::Ref{_mean}, RooAbsReal::Ref{_sigma}} {}
30 RooLandau(const char *name, const char *title, RooAbsReal::Ref _x, RooAbsReal::Ref _mean, RooAbsReal::Ref _sigma);
31 RooLandau(const RooLandau& other, const char* name=nullptr);
32 TObject* clone(const char* newname) const override { return new RooLandau(*this,newname); }
33
34 Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK=true) const override;
35 void generateEvent(Int_t code) override;
36
37 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName = nullptr) const override;
38 double analyticalIntegral(Int_t code, const char *rangeName) const override;
39
40 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
41 std::string
42 buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override;
43
44protected:
45
49
50 double evaluate() const override ;
51 void doEval(RooFit::EvalContext &) const override;
52 inline bool canComputeBatchWithCuda() const override { return true; }
53
54private:
55
56 ClassDefOverride(RooLandau,1) // Landau Distribution PDF
57};
58
59#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:68
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
A class to maintain the context for squashing of RooFit models into code.
Landau distribution p.d.f.
Definition RooLandau.h:24
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition RooLandau.cxx:86
void generateEvent(Int_t code) override
Interface for generation of an event using the algorithm corresponding to the specified code.
TObject * clone(const char *newname) const override
Definition RooLandau.h:32
RooLandau(const char *name, const char *title, RooAbsReal &_x, RooAbsReal &_mean, RooAbsReal &_sigma)
Definition RooLandau.h:28
RooRealProxy x
Definition RooLandau.h:46
RooRealProxy sigma
Definition RooLandau.h:48
void translate(RooFit::Detail::CodeSquashContext &ctx) const override
This function defines a translation for each RooAbsReal based object that can be used to express the ...
Definition RooLandau.cxx:66
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition RooLandau.cxx:59
bool canComputeBatchWithCuda() const override
Definition RooLandau.h:52
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Landau distribution.
Definition RooLandau.cxx:73
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
Definition RooLandau.cxx:81
RooRealProxy mean
Definition RooLandau.h:47
double analyticalIntegral(Int_t code, const char *rangeName) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition RooLandau.cxx:91
std::string buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41