ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
RooExponential.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 /**
18 \file RooExponential.cxx
19 \class RooExponential
20 \ingroup Roofit
21 
22 Exponential p.d.f
23 **/
24 
25 #include "RooFit.h"
26 
27 #include "Riostream.h"
28 #include "Riostream.h"
29 #include <math.h>
30 
31 #include "RooExponential.h"
32 #include "RooRealVar.h"
33 
34 using namespace std;
35 
37 
38 
39 ////////////////////////////////////////////////////////////////////////////////
40 
41 RooExponential::RooExponential(const char *name, const char *title,
42  RooAbsReal& _x, RooAbsReal& _c) :
43  RooAbsPdf(name, title),
44  x("x","Dependent",this,_x),
45  c("c","Exponent",this,_c)
46 {
47 }
48 
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 
53  RooAbsPdf(other, name), x("x",this,other.x), c("c",this,other.c)
54 {
55 }
56 
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 ///cout << "exp(x=" << x << ",c=" << c << ")=" << exp(c*x) << endl ;
60 
62  return exp(c*x);
63 }
64 
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 
68 Int_t RooExponential::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
69 {
70  if (matchArgs(allVars,analVars,x)) return 1 ;
71  return 0 ;
72 }
73 
74 
75 ////////////////////////////////////////////////////////////////////////////////
76 
77 Double_t RooExponential::analyticalIntegral(Int_t code, const char* rangeName) const
78 {
79  switch(code) {
80  case 1:
81  {
82  Double_t ret(0) ;
83  if(c == 0.0) {
84  ret = (x.max(rangeName) - x.min(rangeName));
85  } else {
86  ret = ( exp( c*x.max(rangeName) ) - exp( c*x.min(rangeName) ) )/c;
87  }
88 
89  //cout << "Int_exp_dx(c=" << c << ", xmin=" << x.min(rangeName) << ", xmax=" << x.max(rangeName) << ")=" << ret << endl ;
90  return ret ;
91  }
92  }
93 
94  assert(0) ;
95  return 0 ;
96 }
97 
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Double_t evaluate() const
cout << "exp(x=" << x << ",c=" << c << ")=" << exp(c*x) << endl ;
return c
#define assert(cond)
Definition: unittest.h:542
int Int_t
Definition: RtypesCore.h:41
RooRealProxy x
RooRealProxy c
Double_t x[n]
Definition: legend1.C:17
Exponential p.d.f.
TPaveLabel title(3, 27.1, 15, 28.7,"ROOT Environment and Tools")
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
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
#define name(a, b)
Definition: linkTestLib0.cpp:5
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
ClassImp(RooExponential) RooExponential
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
double exp(double)
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57