Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooCBShape.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id: RooCBShape.h,v 1.11 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_CB_SHAPE
17#define ROO_CB_SHAPE
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21
22class RooRealVar;
23
24class RooCBShape : public RooAbsPdf {
25public:
27 RooCBShape(const char *name, const char *title, RooAbsReal& _m,
28 RooAbsReal& _m0, RooAbsReal& _sigma,
29 RooAbsReal& _alpha, RooAbsReal& _n);
30
31 RooCBShape(const RooCBShape& other, const char *name = nullptr);
32 TObject* clone(const char* newname) const override { return new RooCBShape(*this,newname); }
33
34 Int_t getAnalyticalIntegral( RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr ) const override;
35 double analyticalIntegral(Int_t, const char *rangeName = nullptr) const override;
36
37 // Optimized accept/reject generator support
38 Int_t getMaxVal(const RooArgSet& vars) const override ;
39 double maxVal(Int_t code) const override ;
40
41 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
42 std::string
43 buildCallToAnalyticIntegral(Int_t code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override;
44
45protected:
46
47 double ApproxErf(double arg) const ;
48
54
55 double evaluate() const override;
56 void doEval(RooFit::EvalContext &) const override;
57 inline bool canComputeBatchWithCuda() const override { return true; }
58
59
60private:
61
62 ClassDefOverride(RooCBShape,1) // Crystal Ball lineshape PDF
63};
64
65#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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:55
PDF implementing the Crystal Ball line shape.
Definition RooCBShape.h:24
double ApproxErf(double arg) const
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Crystal ball Shape distribution.
double analyticalIntegral(Int_t, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy n
Definition RooCBShape.h:53
bool canComputeBatchWithCuda() const override
Definition RooCBShape.h:57
RooRealProxy m0
Definition RooCBShape.h:50
RooRealProxy m
Definition RooCBShape.h:49
RooRealProxy sigma
Definition RooCBShape.h:51
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
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 ...
RooRealProxy alpha
Definition RooCBShape.h:52
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise that we know the maximum of self for given (m0,alpha,n,sigma)
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
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.
TObject * clone(const char *newname) const override
Definition RooCBShape.h:32
A class to maintain the context for squashing of RooFit models into code.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41