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 <cmath>
27
28#include "RooBreitWigner.h"
29#include "RooRealVar.h"
30#include "RooHelpers.h"
31#include "RooBatchCompute.h"
32
33
34////////////////////////////////////////////////////////////////////////////////
35
36RooBreitWigner::RooBreitWigner(const char *name, const char *title,
37 RooAbsReal& _x, RooAbsReal& _mean,
39 RooAbsPdf(name,title),
40 x("x","Dependent",this,_x),
41 mean("mean","Mean",this,_mean),
42 width("width","Width",this,_width)
43{
45}
46
47////////////////////////////////////////////////////////////////////////////////
48
50 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
51 width("width",this,other.width)
52{
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
58{
59 double arg= x - mean;
60 return 1. / (arg*arg + 0.25*width*width);
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Compute multiple values of BreitWigner distribution.
66{
67 RooBatchCompute::compute(ctx.config(this), RooBatchCompute::BreitWigner, ctx.output(), {ctx.at(x), ctx.at(mean), ctx.at(width)});
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
72Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
73{
74 if (matchArgs(allVars,analVars,x)) return 1 ;
75 return 0 ;
76}
77
78////////////////////////////////////////////////////////////////////////////////
79
80double RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
81{
82 switch(code) {
83 case 1:
84 {
85 double c = 2./width;
86 return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
87 }
88 }
89
90 assert(0) ;
91 return 0 ;
92}
#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:32
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
Definition RooAbsReal.h:428
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={})
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.