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 "RooFit.h"
32#include "RooRealVar.h"
33#include "RooHelpers.h"
34#include "RooBatchCompute.h"
35
36#include <cmath>
37using namespace std;
38
40
41////////////////////////////////////////////////////////////////////////////////
42/// Construct a Bukin PDF.
43/// \param name The name of the PDF for RooFit's bookkeeping.
44/// \param title The title for e.g. plotting it.
45/// \param _x The variable.
46/// \param _Xp The peak position.
47/// \param _sigp The peak width as FWHM divided by 2*sqrt(2*log(2))=2.35
48/// \param _xi Peak asymmetry. Use values around 0.
49/// \param _rho1 Left tail. Use slightly negative starting values.
50/// \param _rho2 Right tail. Use slightly positive starting values.
51RooBukinPdf::RooBukinPdf(const char *name, const char *title,
52 RooAbsReal& _x, RooAbsReal& _Xp,
53 RooAbsReal& _sigp, RooAbsReal& _xi,
54 RooAbsReal& _rho1, RooAbsReal& _rho2) :
55 // The two addresses refer to our first dependent variable and
56 // parameter, respectively, as declared in the rdl file
57 RooAbsPdf(name, title),
58 x("x","x",this,_x),
59 Xp("Xp","Xp",this,_Xp),
60 sigp("sigp","sigp",this,_sigp),
61 xi("xi","xi",this,_xi),
62 rho1("rho1","rho1",this,_rho1),
63 rho2("rho2","rho2",this,_rho2)
64{
65 RooHelpers::checkRangeOfParameters(this, {&_sigp}, 0.0);
66 RooHelpers::checkRangeOfParameters(this, {&_rho1},-1.0, 0.0);
67 RooHelpers::checkRangeOfParameters(this, {&_rho2}, 0.0, 1.0);
68 RooHelpers::checkRangeOfParameters(this, {&_xi}, -1.0, 1.0);
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// Copy a Bukin PDF.
73RooBukinPdf::RooBukinPdf(const RooBukinPdf& other, const char *name):
74 RooAbsPdf(other,name),
75 x("x",this,other.x),
76 Xp("Xp",this,other.Xp),
77 sigp("sigp",this,other.sigp),
78 xi("xi",this,other.xi),
79 rho1("rho1",this,other.rho1),
80 rho2("rho2",this,other.rho2)
81
82{
83}
84
85////////////////////////////////////////////////////////////////////////////////
86/// Implementation
87
89{
90 const double consts = 2*sqrt(2*log(2.0));
91 double r1=0,r2=0,r3=0,r4=0,r5=0,hp=0;
92 double x1 = 0,x2 = 0;
93 double fit_result = 0;
94
95 hp=sigp*consts;
96 r3=log(2.);
97 r4=sqrt(xi*xi+1);
98 r1=xi/r4;
99
100 if(fabs(xi) > exp(-6.)){
101 r5=xi/log(r4+xi);
102 }
103 else
104 r5=1;
105
106 x1 = Xp + (hp / 2) * (r1-1);
107 x2 = Xp + (hp / 2) * (r1+1);
108
109 //--- Left Side
110 if(x < x1){
111 r2=rho1*(x-x1)*(x-x1)/(Xp-x1)/(Xp-x1)-r3 + 4 * r3 * (x-x1)/hp * r5 * r4/(r4-xi)/(r4-xi);
112 }
113
114
115 //--- Center
116 else if(x < x2) {
117 if(fabs(xi) > exp(-6.)) {
118 r2=log(1 + 4 * xi * r4 * (x-Xp)/hp)/log(1+2*xi*(xi-r4));
119 r2=-r3*r2*r2;
120 }
121 else{
122 r2=-4*r3*(x-Xp)*(x-Xp)/hp/hp;
123 }
124 }
125
126
127 //--- Right Side
128 else {
129 r2=rho2*(x-x2)*(x-x2)/(Xp-x2)/(Xp-x2)-r3 - 4 * r3 * (x-x2)/hp * r5 * r4/(r4+xi)/(r4+xi);
130 }
131
132 if(fabs(r2) > 100){
133 fit_result = 0;
134 }
135 else{
136 //---- Normalize the result
137 fit_result = exp(r2);
138 }
139
140 return fit_result;
141}
142
143////////////////////////////////////////////////////////////////////////////////
144/// Compute multiple values of Bukin distribution.
145void RooBukinPdf::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
146{
148 dispatch->compute(stream, RooBatchCompute::Bukin, output, nEvents,
149 {dataMap.at(x), dataMap.at(Xp), dataMap.at(sigp), dataMap.at(xi), dataMap.at(rho1), dataMap.at(rho2)});
150}
static const double x2[5]
static const double x1[5]
#define ClassImp(name)
Definition Rtypes.h:364
char name[80]
Definition TGX11.cxx:110
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
The RooBukinPdf implements the NovosibirskA function.
Definition RooBukinPdf.h:29
RooRealProxy rho2
Definition RooBukinPdf.h:49
RooRealProxy sigp
Definition RooBukinPdf.h:46
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const
Compute multiple values of Bukin distribution.
RooRealProxy rho1
Definition RooBukinPdf.h:48
RooRealProxy x
Definition RooBukinPdf.h:44
double evaluate() const
Implementation.
RooRealProxy xi
Definition RooBukinPdf.h:47
RooRealProxy Xp
Definition RooBukinPdf.h:45
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
Double_t x[n]
Definition legend1.C:17
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
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.
static void output(int code)
Definition gifencode.c:226