Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 <cmath>
27
28#include "RooBreitWigner.h"
29#include "RooRealVar.h"
30#include "RooBatchCompute.h"
31
32
33////////////////////////////////////////////////////////////////////////////////
34
35RooBreitWigner::RooBreitWigner(const char *name, const char *title,
36 RooAbsReal& _x, RooAbsReal& _mean,
38 RooAbsPdf(name,title),
39 x("x","Dependent",this,_x),
40 mean("mean","Mean",this,_mean),
41 width("width","Width",this,_width)
42{
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
48 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
49 width("width",this,other.width)
50{
51}
52
53////////////////////////////////////////////////////////////////////////////////
54
56{
57 double arg= x - mean;
58 return 1. / (arg*arg + 0.25*width*width);
59}
60
61////////////////////////////////////////////////////////////////////////////////
62/// Compute multiple values of BreitWigner distribution.
64{
65 RooBatchCompute::compute(ctx.config(this), RooBatchCompute::BreitWigner, ctx.output(), {ctx.at(x), ctx.at(mean), ctx.at(width)});
66}
67
68////////////////////////////////////////////////////////////////////////////////
69
70Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
71{
72 if (matchArgs(allVars,analVars,x)) return 1 ;
73 return 0 ;
74}
75
76////////////////////////////////////////////////////////////////////////////////
77
78double RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
79{
80 switch(code) {
81 case 1:
82 {
83 double c = 2./width;
84 return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
85 }
86 }
87
88 assert(0) ;
89 return 0 ;
90}
#define c(i)
Definition RSha256.hxx:101
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t width
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
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
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.
RooRealProxy mean
void doEval(RooFit::EvalContext &) const override
Compute multiple values of BreitWigner distribution.
RooRealProxy x
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
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
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})