ROOT  6.06/09
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 //////////////////////////////////////////////////////////////////////////////
19 //
20 // BEGIN_HTML
21 // Class RooBreitWigner is a RooAbsPdf implementation
22 // that models a non-relativistic Breit-Wigner shape
23 // END_HTML
24 //
25 
26 
27 #include "RooFit.h"
28 
29 #include "Riostream.h"
30 #include "Riostream.h"
31 #include <math.h>
32 
33 #include "RooBreitWigner.h"
34 #include "RooAbsReal.h"
35 #include "RooRealVar.h"
36 // #include "RooFitTools/RooRandom.h"
37 
38 using namespace std;
39 
41 
42 
43 ////////////////////////////////////////////////////////////////////////////////
44 
45 RooBreitWigner::RooBreitWigner(const char *name, const char *title,
46  RooAbsReal& _x, RooAbsReal& _mean,
47  RooAbsReal& _width) :
48  RooAbsPdf(name,title),
49  x("x","Dependent",this,_x),
50  mean("mean","Mean",this,_mean),
51  width("width","Width",this,_width)
52 {
53 }
54 
55 
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 
59 RooBreitWigner::RooBreitWigner(const RooBreitWigner& other, const char* name) :
60  RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
61  width("width",this,other.width)
62 {
63 }
64 
65 
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 
70 {
71  Double_t arg= x - mean;
72  return 1. / (arg*arg + 0.25*width*width);
73 }
74 
75 
76 
77 ////////////////////////////////////////////////////////////////////////////////
78 
79 Int_t RooBreitWigner::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
80 {
81  if (matchArgs(allVars,analVars,x)) return 1 ;
82  return 0 ;
83 }
84 
85 
86 
87 ////////////////////////////////////////////////////////////////////////////////
88 
89 Double_t RooBreitWigner::analyticalIntegral(Int_t code, const char* rangeName) const
90 {
91  switch(code) {
92  case 1:
93  {
94  Double_t c = 2./width;
95  return c*(atan(c*(x.max(rangeName)-mean)) - atan(c*(x.min(rangeName)-mean)));
96  }
97  }
98 
99  assert(0) ;
100  return 0 ;
101 }
102 
#define assert(cond)
Definition: unittest.h:542
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Double_t x[n]
Definition: legend1.C:17
RooRealProxy mean
RooRealProxy width
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:824
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Double_t evaluate() const
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
ClassImp(RooBreitWigner) RooBreitWigner
double atan(double)
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooRealProxy x
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57