Logo ROOT   6.16/01
Reference Guide
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 "RooFit.h"
26
27#include "Riostream.h"
28#include "Riostream.h"
29#include <math.h>
30
31#include "RooBreitWigner.h"
32#include "RooAbsReal.h"
33#include "RooRealVar.h"
34// #include "RooFitTools/RooRandom.h"
35
36using namespace std;
37
39
40////////////////////////////////////////////////////////////////////////////////
41
42RooBreitWigner::RooBreitWigner(const char *name, const char *title,
43 RooAbsReal& _x, RooAbsReal& _mean,
44 RooAbsReal& _width) :
45 RooAbsPdf(name,title),
46 x("x","Dependent",this,_x),
47 mean("mean","Mean",this,_mean),
48 width("width","Width",this,_width)
49{
50}
51
52////////////////////////////////////////////////////////////////////////////////
53
55 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
56 width("width",this,other.width)
57{
58}
59
60////////////////////////////////////////////////////////////////////////////////
61
63{
64 Double_t arg= x - mean;
65 return 1. / (arg*arg + 0.25*width*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_t RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
79{
80 switch(code) {
81 case 1:
82 {
83 Double_t 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
int Int_t
Definition: RtypesCore.h:41
double Double_t
Definition: RtypesCore.h:55
#define ClassImp(name)
Definition: Rtypes.h:363
include TDocParser_001 C image html pict1_TDocParser_001 png width
Definition: TDocParser.cxx:121
double atan(double)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
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
Class RooBreitWigner is a RooAbsPdf implementation that models a non-relativistic Breit-Wigner shape.
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
RooRealProxy width
RooRealProxy mean
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy x
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
Double_t x[n]
Definition: legend1.C:17
STL namespace.