Logo ROOT   6.14/05
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 
21 Class RooBreitWigner is a RooAbsPdf implementation
22 that 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 
36 using namespace std;
37 
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
42 RooBreitWigner::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 
70 Int_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 
78 Double_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 }
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
image html pict1_TGaxis_012 png width
Define new text attributes for the label number "labNum".
Definition: TGaxis.cxx:2551
Double_t evaluate() const
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
int Int_t
Definition: RtypesCore.h:41
STL namespace.
Double_t x[n]
Definition: legend1.C:17
RooRealProxy mean
RooRealProxy width
Class RooBreitWigner is a RooAbsPdf implementation that models a non-relativistic Breit-Wigner shape...
#define ClassImp(name)
Definition: Rtypes.h:359
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
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
double atan(double)
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
RooRealProxy x
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
#define c(i)
Definition: RSha256.hxx:101
char name[80]
Definition: TGX11.cxx:109