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 "RooFit.h"
26
27#include "Riostream.h"
28#include <math.h>
29
30#include "RooBreitWigner.h"
31#include "RooAbsReal.h"
32#include "RooRealVar.h"
33#include "RooBatchCompute.h"
34
35using namespace std;
36
38
39////////////////////////////////////////////////////////////////////////////////
40
41RooBreitWigner::RooBreitWigner(const char *name, const char *title,
42 RooAbsReal& _x, RooAbsReal& _mean,
43 RooAbsReal& _width) :
44 RooAbsPdf(name,title),
45 x("x","Dependent",this,_x),
46 mean("mean","Mean",this,_mean),
47 width("width","Width",this,_width)
48{
49}
50
51////////////////////////////////////////////////////////////////////////////////
52
54 RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
55 width("width",this,other.width)
56{
57}
58
59////////////////////////////////////////////////////////////////////////////////
60
62{
63 Double_t arg= x - mean;
64 return 1. / (arg*arg + 0.25*width*width);
65}
66
67////////////////////////////////////////////////////////////////////////////////
68/// Compute multiple values of BreitWigner distribution.
70 return RooBatchCompute::dispatch->computeBreitWigner(this, evalData, x->getValues(evalData, normSet), mean->getValues(evalData, normSet), width->getValues(evalData, normSet));
71}
72
73////////////////////////////////////////////////////////////////////////////////
74
75Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
76{
77 if (matchArgs(allVars,analVars,x)) return 1 ;
78 return 0 ;
79}
80
81////////////////////////////////////////////////////////////////////////////////
82
83Double_t RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
84{
85 switch(code) {
86 case 1:
87 {
88 Double_t c = 2./width;
89 return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
90 }
91 }
92
93 assert(0) ;
94 return 0 ;
95}
#define c(i)
Definition RSha256.hxx:101
#define ClassImp(name)
Definition Rtypes.h:364
include TDocParser_001 C image html pict1_TDocParser_001 png width
char name[80]
Definition TGX11.cxx:110
double atan(double)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
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:29
virtual RooSpan< double > computeBreitWigner(const RooAbsReal *, RunContext &, RooSpan< const double > x, RooSpan< const double > mean, RooSpan< const double > width)=0
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.
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const
Compute multiple values of BreitWigner distribution.
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
A simple container to hold a batch of data values.
Definition RooSpan.h:34
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
Double_t x[n]
Definition legend1.C:17
R__EXTERN RooBatchComputeInterface * dispatch
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31