Logo ROOT   6.18/05
Reference Guide
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 "RooFit.h"
29
30#include "Riostream.h"
31#include "Riostream.h"
32#include <math.h>
33
34#include "RooArgusBG.h"
35#include "RooRealVar.h"
36#include "RooRealConstant.h"
37#include "RooMath.h"
38#include "TMath.h"
39
40#include "TError.h"
41
42using namespace std;
43
45
46////////////////////////////////////////////////////////////////////////////////
47/// Constructor.
48
49RooArgusBG::RooArgusBG(const char *name, const char *title,
50 RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _c) :
51 RooAbsPdf(name, title),
52 m("m","Mass",this,_m),
53 m0("m0","Resonance mass",this,_m0),
54 c("c","Slope parameter",this,_c),
55 p("p","Power",this,(RooRealVar&)RooRealConstant::value(0.5))
56{
57}
58
59////////////////////////////////////////////////////////////////////////////////
60/// Constructor.
61
62RooArgusBG::RooArgusBG(const char *name, const char *title,
63 RooAbsReal& _m, RooAbsReal& _m0, RooAbsReal& _c, RooAbsReal& _p) :
64 RooAbsPdf(name, title),
65 m("m","Mass",this,_m),
66 m0("m0","Resonance mass",this,_m0),
67 c("c","Slope parameter",this,_c),
68 p("p","Power",this,_p)
69{
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// Constructor.
74
75RooArgusBG::RooArgusBG(const RooArgusBG& other, const char* name) :
76 RooAbsPdf(other,name),
77 m("m",this,other.m),
78 m0("m0",this,other.m0),
79 c("c",this,other.c),
80 p("p",this,other.p)
81{
82}
83
84////////////////////////////////////////////////////////////////////////////////
85
87 Double_t t= m/m0;
88 if(t >= 1) return 0;
89
90 Double_t u= 1 - t*t;
91 //cout << "c = " << c << " result = " << m*TMath::Power(u,p)*exp(c*u) << endl ;
92 return m*TMath::Power(u,p)*exp(c*u) ;
93}
94
95////////////////////////////////////////////////////////////////////////////////
96
97Int_t RooArgusBG::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
98{
99 if (p.arg().isConstant()) {
100 // We can integrate over m if power = 0.5
101 if (matchArgs(allVars,analVars,m) && p == 0.5) return 1;
102 }
103 return 0;
104
105}
106
107////////////////////////////////////////////////////////////////////////////////
108
109Double_t RooArgusBG::analyticalIntegral(Int_t code, const char* rangeName) const
110{
111 R__ASSERT(code==1);
112 // Formula for integration over m when p=0.5
113 static const Double_t pi = atan2(0.0,-1.0);
114 Double_t min = (m.min(rangeName) < m0) ? m.min(rangeName) : m0;
115 Double_t max = (m.max(rangeName) < m0) ? m.max(rangeName) : m0;
116 Double_t f1 = (1.-TMath::Power(min/m0,2));
117 Double_t f2 = (1.-TMath::Power(max/m0,2));
118 Double_t aLow, aHigh ;
119 if ( c < 0. ) {
120 aLow = -0.5*m0*m0*(exp(c*f1)*sqrt(f1)/c + 0.5/TMath::Power(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f1)));
121 aHigh = -0.5*m0*m0*(exp(c*f2)*sqrt(f2)/c + 0.5/TMath::Power(-c,1.5)*sqrt(pi)*RooMath::erf(sqrt(-c*f2)));
122 } else if ( c == 0. ) {
123 aLow = -m0*m0/3.*f1*sqrt(f1);
124 aHigh = -m0*m0/3.*f1*sqrt(f2);
125 } else {
126 aLow = 0.5*m0*m0*exp(c*f1)/(c*sqrt(c)) * (0.5*sqrt(pi)*(RooMath::faddeeva(sqrt(c*f1))).imag() - sqrt(c*f1));
127 aHigh = 0.5*m0*m0*exp(c*f2)/(c*sqrt(c)) * (0.5*sqrt(pi)*(RooMath::faddeeva(sqrt(c*f2))).imag() - sqrt(c*f2));
128 }
129 Double_t area = aHigh - aLow;
130 //cout << "c = " << c << "aHigh = " << aHigh << " aLow = " << aLow << " area = " << area << endl ;
131 return area;
132
133}
#define c(i)
Definition: RSha256.hxx:101
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:365
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
double atan2(double, double)
double sqrt(double)
double exp(double)
Bool_t isConstant() const
Definition: RooAbsArg.h:311
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Bool_t 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:28
RooArgusBG is a RooAbsPdf implementation describing the ARGUS background shape.
Definition: RooArgusBG.h:25
RooRealProxy m
Definition: RooArgusBG.h:40
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Definition: RooArgusBG.cxx:86
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Definition: RooArgusBG.cxx:97
RooRealProxy c
Definition: RooArgusBG.h:42
RooRealProxy p
Definition: RooArgusBG.h:43
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Definition: RooArgusBG.cxx:109
RooRealProxy m0
Definition: RooArgusBG.h:41
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:580
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition: RooMath.cxx:542
RooRealConstant provides static functions to create and keep track of RooRealVar constants.
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
const RooAbsReal & arg() const
Definition: RooRealProxy.h:43
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
TF1 * f1
Definition: legend1.C:11
static constexpr double pi
LongDouble_t Power(LongDouble_t x, LongDouble_t y)
Definition: TMath.h:723
auto * m
Definition: textangle.C:8