Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBukinPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * RW, Ruddick William UC Colorado wor@slac.stanford.edu *
7 * *
8 * Copyright (c) 2000-2005, Regents of the University of California *
9 * and Stanford University. All rights reserved. *
10 * *
11 * Redistribution and use in source and binary forms, *
12 * with or without modification, are permitted according to the terms *
13 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
14 *****************************************************************************/
15
16/** \class RooBukinPdf
17 \ingroup Roofit
18
19The RooBukinPdf implements the NovosibirskA function. For the parameters, see
20RooBukinPdf().
21
22Credits:
23May 26, 2003.
24A.Bukin, Budker INP, Novosibirsk
25
26\image html RooBukin.png
27http://www.slac.stanford.edu/BFROOT/www/Organization/CollabMtgs/2003/detJuly2003/Tues3a/bukin.ps
28**/
29
30#include "RooBukinPdf.h"
31#include "RooRealVar.h"
32#include "RooHelpers.h"
33#include "RooBatchCompute.h"
34
35#include <cmath>
36
38
39////////////////////////////////////////////////////////////////////////////////
40/// Construct a Bukin PDF.
41/// \param name The name of the PDF for RooFit's bookkeeping.
42/// \param title The title for e.g. plotting it.
43/// \param _x The variable.
44/// \param _Xp The peak position.
45/// \param _sigp The peak width as FWHM divided by 2*sqrt(2*log(2))=2.35
46/// \param _xi Peak asymmetry. Use values around 0.
47/// \param _rho1 Left tail. Use slightly negative starting values.
48/// \param _rho2 Right tail. Use slightly positive starting values.
49RooBukinPdf::RooBukinPdf(const char *name, const char *title,
50 RooAbsReal& _x, RooAbsReal& _Xp,
51 RooAbsReal& _sigp, RooAbsReal& _xi,
52 RooAbsReal& _rho1, RooAbsReal& _rho2) :
53 // The two addresses refer to our first dependent variable and
54 // parameter, respectively, as declared in the rdl file
55 RooAbsPdf(name, title),
56 x("x","x",this,_x),
57 Xp("Xp","Xp",this,_Xp),
58 sigp("sigp","sigp",this,_sigp),
59 xi("xi","xi",this,_xi),
60 rho1("rho1","rho1",this,_rho1),
61 rho2("rho2","rho2",this,_rho2)
62{
63 RooHelpers::checkRangeOfParameters(this, {&_sigp}, 0.0);
64 RooHelpers::checkRangeOfParameters(this, {&_rho1},-1.0, 0.0);
65 RooHelpers::checkRangeOfParameters(this, {&_rho2}, 0.0, 1.0);
66 RooHelpers::checkRangeOfParameters(this, {&_xi}, -1.0, 1.0);
67}
68
69////////////////////////////////////////////////////////////////////////////////
70/// Copy a Bukin PDF.
71RooBukinPdf::RooBukinPdf(const RooBukinPdf& other, const char *name):
72 RooAbsPdf(other,name),
73 x("x",this,other.x),
74 Xp("Xp",this,other.Xp),
75 sigp("sigp",this,other.sigp),
76 xi("xi",this,other.xi),
77 rho1("rho1",this,other.rho1),
78 rho2("rho2",this,other.rho2)
79
80{
81}
82
83////////////////////////////////////////////////////////////////////////////////
84/// Implementation
85
87{
88 const double consts = 2*sqrt(2*log(2.0));
89 double r1 = 0;
90 double r2 = 0;
91 double r3 = 0;
92 double r4 = 0;
93 double r5 = 0;
94 double hp = 0;
95 double x1 = 0;
96 double x2 = 0;
97 double fit_result = 0;
98
99 hp=sigp*consts;
100 r3=log(2.);
101 r4=sqrt(xi*xi+1);
102 r1=xi/r4;
103
104 if(std::abs(xi) > exp(-6.)){
105 r5=xi/log(r4+xi);
106 }
107 else
108 r5=1;
109
110 x1 = Xp + (hp / 2) * (r1-1);
111 x2 = Xp + (hp / 2) * (r1+1);
112
113 //--- Left Side
114 if(x < x1){
115 r2=rho1*(x-x1)*(x-x1)/(Xp-x1)/(Xp-x1)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/(r4-xi)/(r4-xi);
116 }
117
118
119 //--- Center
120 else if(x < x2) {
121 if(std::abs(xi) > exp(-6.)) {
122 r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
123 r2=-r3*r2*r2;
124 }
125 else{
126 r2=-4*r3*(x-Xp)*(x-Xp)/hp/hp;
127 }
128 }
129
130
131 //--- Right Side
132 else {
133 r2=rho2*(x-x2)*(x-x2)/(Xp-x2)/(Xp-x2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/(r4+xi)/(r4+xi);
134 }
135
136 if(std::abs(r2) > 100){
137 fit_result = 0;
138 }
139 else{
140 //---- Normalize the result
141 fit_result = exp(r2);
142 }
143
144 return fit_result;
145}
146
147////////////////////////////////////////////////////////////////////////////////
148/// Compute multiple values of Bukin distribution.
150{
152 {ctx.at(x), ctx.at(Xp), ctx.at(sigp), ctx.at(xi), ctx.at(rho1), ctx.at(rho2)});
153}
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
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
The RooBukinPdf implements the NovosibirskA function.
Definition RooBukinPdf.h:29
RooRealProxy rho2
Definition RooBukinPdf.h:48
RooRealProxy sigp
Definition RooBukinPdf.h:45
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Bukin distribution.
RooRealProxy rho1
Definition RooBukinPdf.h:47
RooRealProxy x
Definition RooBukinPdf.h:43
RooRealProxy xi
Definition RooBukinPdf.h:46
double evaluate() const override
Implementation.
RooRealProxy Xp
Definition RooBukinPdf.h:44
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string const &extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.