Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBreitWigner.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * AS, Abi Soffer, Colorado State University, abi@slac.stanford.edu *
7 * TS, Thomas Schietinger, SLAC, schieti@slac.stanford.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * Colorado State University *
11 * and Stanford University. All rights reserved. *
12 * *
13 * Redistribution and use in source and binary forms, *
14 * with or without modification, are permitted according to the terms *
15 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
16 *****************************************************************************/
17
18/** \class RooBreitWigner
19 \ingroup Roofit
20
21Class RooBreitWigner is a RooAbsPdf implementation
22that models a non-relativistic Breit-Wigner shape
23**/
24
25#include "Riostream.h"
26#include <math.h>
27
28#include "RooBreitWigner.h"
29#include "RooRealVar.h"
30#include "RooBatchCompute.h"
31
32using namespace std;
33
35
36////////////////////////////////////////////////////////////////////////////////
37
38RooBreitWigner::RooBreitWigner(const char *name, const char *title,
39 RooAbsReal& _x, RooAbsReal& _mean,
40 RooAbsReal& _width) :
41 RooAbsPdf(name,title),
42 x("x","Dependent",this,_x),
43 mean("mean","Mean",this,_mean),
44 width("width","Width",this,_width)
45{
46}
47
48////////////////////////////////////////////////////////////////////////////////
49
51 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
52 width("width",this,other.width)
53{
54}
55
56////////////////////////////////////////////////////////////////////////////////
57
59{
60 double arg= x - mean;
61 return 1. / (arg*arg + 0.25*width*width);
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// Compute multiple values of BreitWigner distribution.
66void RooBreitWigner::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
67{
69 dispatch->compute(stream, RooBatchCompute::BreitWigner, output, nEvents, {dataMap.at(x), dataMap.at(mean), dataMap.at(width)});
70}
71
72////////////////////////////////////////////////////////////////////////////////
73
74Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
75{
76 if (matchArgs(allVars,analVars,x)) return 1 ;
77 return 0 ;
78}
79
80////////////////////////////////////////////////////////////////////////////////
81
82double RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
83{
84 switch(code) {
85 case 1:
86 {
87 double c = 2./width;
88 return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
89 }
90 }
91
92 assert(0) ;
93 return 0 ;
94}
#define c(i)
Definition RSha256.hxx:101
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t width
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:62
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:55
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, ArgVector &)=0
Class RooBreitWigner is a RooAbsPdf implementation that models a non-relativistic Breit-Wigner shape.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy width
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void computeBatch(cudaStream_t *, double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of BreitWigner distribution.
RooRealProxy mean
RooRealProxy x
RooSpan< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:21
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
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...
static void output()