Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
RooCrystalBall.h
Go to the documentation of this file.
1// Author: Jonas Rembser, CERN 02/2021
2
3#ifndef RooFit_RooFit_RooCrystalBall_h
4#define RooFit_RooFit_RooCrystalBall_h
5
6#include "RooAbsPdf.h"
7#include "RooRealProxy.h"
8
9#include <memory>
10
12public:
14
15 RooCrystalBall(const char *name, const char *title, RooAbsReal &x, RooAbsReal &x0, RooAbsReal &sigmaL,
17 RooCrystalBall(const char *name, const char *title, RooAbsReal &x, RooAbsReal &x0, RooAbsReal &sigmaLR,
19 RooCrystalBall(const char *name, const char *title, RooAbsReal &x, RooAbsReal &x0, RooAbsReal &sigmaLR,
20 RooAbsReal &alpha, RooAbsReal &n, bool doubleSided = false);
21
22 RooCrystalBall(const RooCrystalBall &other, const char *name = nullptr);
23 TObject *clone(const char *newname) const override { return new RooCrystalBall(*this, newname); }
24
25 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName = nullptr) const override;
26 double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override;
27
28 // Optimized accept/reject generator support
29 Int_t getMaxVal(const RooArgSet &vars) const override;
30 double maxVal(Int_t code) const override;
31
32 // Getters for non-optional parameters
33 RooAbsReal const &x() const { return *x_; }
34 RooAbsReal const &x0() const { return *x0_; }
35 RooAbsReal const &sigmaL() const { return *sigmaL_; }
36 RooAbsReal const &sigmaR() const { return *sigmaR_; }
37 RooAbsReal const &alphaL() const { return *alphaL_; }
38 RooAbsReal const &nL() const { return *nL_; }
39
40 // Getters for optional parameter: return nullptr if parameter is not set
41 RooAbsReal const *alphaR() const { return alphaR_ ? &**alphaR_ : nullptr; }
42 RooAbsReal const *nR() const { return nR_ ? &**nR_ : nullptr; }
43
44 // Convenience functions to check if optional parameters are set
45 bool hasAlphaR() const { return alphaR_ != nullptr; }
46 bool hasNR() const { return nR_ != nullptr; }
47
48protected:
49 double evaluate() const override;
50
51private:
58
59 // optional parameters
60 std::unique_ptr<RooRealProxy> alphaR_ = nullptr;
61 std::unique_ptr<RooRealProxy> nR_ = nullptr;
62
64};
65
66#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
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
PDF implementing the generalized Asymmetrical Double-Sided Crystall Ball line shape.
TObject * clone(const char *newname) const override
std::unique_ptr< RooRealProxy > nR_
bool hasAlphaR() const
bool hasNR() const
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooAbsReal const * nR() const
RooAbsReal const * alphaR() const
RooRealProxy sigmaR_
RooRealProxy nL_
RooAbsReal const & sigmaR() const
RooRealProxy x_
std::unique_ptr< RooRealProxy > alphaR_
RooAbsReal const & nL() const
RooRealProxy x0_
RooAbsReal const & x0() const
RooAbsReal const & x() const
RooRealProxy alphaL_
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooAbsReal const & sigmaL() const
RooRealProxy sigmaL_
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise that we know the maximum of self for given (m0,alpha,n,sigma).
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooAbsReal const & alphaL() const
Mother of all ROOT objects.
Definition TObject.h:41
const Int_t n
Definition legend1.C:16