Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooArgusBG.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
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
17/** \class RooArgusBG
18 \ingroup Roofit
19
20RooArgusBG is a RooAbsPdf implementation describing the ARGUS background shape.
21\f[
22 \mathrm{Argus}(m, m_0, c, p) = \mathcal{N} \cdot m \cdot \left[ 1 - \left( \frac{m}{m_0} \right)^2 \right]^p
23 \cdot \exp\left[ c \cdot \left(1 - \left(\frac{m}{m_0}\right)^2 \right) \right]
24\f]
25\image html RooArgusBG.png
26*/
27
28#include "RooArgusBG.h"
29#include "RooRealVar.h"
30#include "RooRealConstant.h"
31#include "RooMath.h"
32#include "RooBatchCompute.h"
33
34#include "TMath.h"
35
36#include <cmath>
37
39
40////////////////////////////////////////////////////////////////////////////////
41/// Constructor.
42
43RooArgusBG::RooArgusBG(const char *name, const char *title,
45 RooAbsPdf(name, title),
46 m("m","Mass",this,_m),
47 m0("m0","Resonance mass",this,_m0),
48 c("c","Slope parameter",this,_c),
49 p("p","Power",this,_p)
50{
51}
52
53////////////////////////////////////////////////////////////////////////////////
54/// Constructor.
55
56RooArgusBG::RooArgusBG(const RooArgusBG& other, const char* name) :
57 RooAbsPdf(other,name),
58 m("m",this,other.m),
59 m0("m0",this,other.m0),
60 c("c",this,other.c),
61 p("p",this,other.p)
62{
63}
64
65////////////////////////////////////////////////////////////////////////////////
66
67double RooArgusBG::evaluate() const {
68 double t= m/m0;
69 if(t >= 1) return 0;
70
71 double u= 1 - t*t;
72 return m*std::pow(u,p)*exp(c*u) ;
73}
74
75////////////////////////////////////////////////////////////////////////////////
76
78{
80 {ctx.at(m), ctx.at(m0), ctx.at(c), ctx.at(p)});
81}
82
83////////////////////////////////////////////////////////////////////////////////
84
85Int_t RooArgusBG::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
86{
87 if (p.arg().isConstant()) {
88 // We can integrate over m if power = 0.5
89 if (matchArgs(allVars,analVars,m) && p == 0.5) return 1;
90 }
91 return 0;
92
93}
94
95////////////////////////////////////////////////////////////////////////////////
96
97double RooArgusBG::analyticalIntegral(Int_t code, const char* rangeName) const
98{
99 R__ASSERT(code==1);
100 // Formula for integration over m when p=0.5
101 static const double pi = atan2(0.0,-1.0);
102 double min = (m.min(rangeName) < m0) ? m.min(rangeName) : m0;
103 double max = (m.max(rangeName) < m0) ? m.max(rangeName) : m0;
104 double f1 = (1.-std::pow(min/m0,2));
105 double f2 = (1.-std::pow(max/m0,2));
106 double aLow;
107 double aHigh;
108 if ( c < 0. ) {
109 aLow = -0.5*m0*m0*(exp(c*f1)*sqrt(f1)/c + 0.5/std::pow(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f1)));
110 aHigh = -0.5*m0*m0*(exp(c*f2)*sqrt(f2)/c + 0.5/std::pow(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f2)));
111 } else if ( c == 0. ) {
112 aLow = -m0*m0/3.*f1*sqrt(f1);
113 aHigh = -m0*m0/3.*f1*sqrt(f2);
114 } else {
115 aLow = 0.5*m0*m0*exp(c*f1)/(c*sqrt(c)) * (0.5*sqrt(pi)*(RooMath::faddeeva(sqrt(c*f1))).imag() - sqrt(c*f1));
116 aHigh = 0.5*m0*m0*exp(c*f2)/(c*sqrt(c)) * (0.5*sqrt(pi)*(RooMath::faddeeva(sqrt(c*f2))).imag() - sqrt(c*f2));
117 }
118 double area = aHigh - aLow;
119 //cout << "c = " << c << "aHigh = " << aHigh << " aLow = " << aLow << " area = " << area << endl ;
120 return area;
121
122}
#define c(i)
Definition RSha256.hxx:101
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
winID h TVirtualViewer3D TVirtualGLPainter p
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
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgusBG is a RooAbsPdf implementation describing the ARGUS background shape.
Definition RooArgusBG.h:22
RooRealProxy m
Definition RooArgusBG.h:42
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy c
Definition RooArgusBG.h:44
RooRealProxy p
Definition RooArgusBG.h:45
RooRealProxy m0
Definition RooArgusBG.h:43
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition RooMath.cxx:59
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition RooMath.cxx:30
TF1 * f1
Definition legend1.C:11
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
TMarker m
Definition textangle.C:8