ROOT  6.06/09
Reference Guide
RooBifurGauss.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * Abi Soffer, Colorado State University, abi@slac.stanford.edu *
7  * *
8  * Copyright (c) 2000-2005, Regents of the University of California, *
9  * Colorado State University *
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 #include "RooFit.h"
17 
18 //////////////////////////////////////////////////////////////////////////////
19 //
20 // BEGIN_HTML
21 // Bifurcated Gaussian p.d.f with different widths on left and right
22 // side of maximum value
23 // END_HTML
24 //
25 
26 
27 #include "Riostream.h"
28 #include "TMath.h"
29 #include <math.h>
30 
31 #include "RooBifurGauss.h"
32 #include "RooAbsReal.h"
33 #include "RooMath.h"
34 
35 using namespace std;
36 
38 
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
42 RooBifurGauss::RooBifurGauss(const char *name, const char *title,
43  RooAbsReal& _x, RooAbsReal& _mean,
44  RooAbsReal& _sigmaL, RooAbsReal& _sigmaR) :
45  RooAbsPdf(name, title),
46  x ("x" , "Dependent" , this, _x),
47  mean ("mean" , "Mean" , this, _mean),
48  sigmaL("sigmaL", "Left Sigma" , this, _sigmaL),
49  sigmaR("sigmaR", "Right Sigma", this, _sigmaR)
50 
51 {
52 }
53 
54 
55 ////////////////////////////////////////////////////////////////////////////////
56 
57 RooBifurGauss::RooBifurGauss(const RooBifurGauss& other, const char* name) :
58  RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
59  sigmaL("sigmaL",this,other.sigmaL), sigmaR("sigmaR", this, other.sigmaR)
60 {
61 }
62 
63 
64 ////////////////////////////////////////////////////////////////////////////////
65 
67  Double_t arg = x - mean;
68 
69  Double_t coef(0.0);
70 
71  if (arg < 0.0){
72  if (TMath::Abs(sigmaL) > 1e-30) {
73  coef = -0.5/(sigmaL*sigmaL);
74  }
75  } else {
76  if (TMath::Abs(sigmaR) > 1e-30) {
77  coef = -0.5/(sigmaR*sigmaR);
78  }
79  }
80 
81  return exp(coef*arg*arg);
82 }
83 
84 
85 ////////////////////////////////////////////////////////////////////////////////
86 
87 Int_t RooBifurGauss::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
88 {
89  if (matchArgs(allVars,analVars,x)) return 1 ;
90  return 0 ;
91 }
92 
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 
96 Double_t RooBifurGauss::analyticalIntegral(Int_t code, const char* rangeName) const
97 {
98  switch(code) {
99  case 1:
100  {
101  static Double_t root2 = sqrt(2.) ;
102  static Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
103 
104 // Double_t coefL(0.0), coefR(0.0);
105 // if (TMath::Abs(sigmaL) > 1e-30) {
106 // coefL = -0.5/(sigmaL*sigmaL);
107 // }
108 
109 // if (TMath::Abs(sigmaR) > 1e-30) {
110 // coefR = -0.5/(sigmaR*sigmaR);
111 // }
112 
113  Double_t xscaleL = root2*sigmaL;
114  Double_t xscaleR = root2*sigmaR;
115 
116  Double_t integral = 0.0;
117  if(x.max(rangeName) < mean)
118  {
119  integral = sigmaL * ( RooMath::erf((x.max(rangeName) - mean)/xscaleL) - RooMath::erf((x.min(rangeName) - mean)/xscaleL) );
120  }
121  else if (x.min(rangeName) > mean)
122  {
123  integral = sigmaR * ( RooMath::erf((x.max(rangeName) - mean)/xscaleR) - RooMath::erf((x.min(rangeName) - mean)/xscaleR) );
124  }
125  else
126  {
127  integral = sigmaR*RooMath::erf((x.max(rangeName) - mean)/xscaleR) - sigmaL*RooMath::erf((x.min(rangeName) - mean)/xscaleL);
128  }
129  // return rootPiBy2*(sigmaR*RooMath::erf((x.max(rangeName) - mean)/xscaleR) -
130  // sigmaL*RooMath::erf((x.min(rangeName) - mean)/xscaleL));
131  return integral*rootPiBy2;
132  }
133  }
134 
135  assert(0) ;
136  return 0 ; // to prevent compiler warnings
137 }
RooRealProxy x
Definition: RooBifurGauss.h:40
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
#define assert(cond)
Definition: unittest.h:542
ClassImp(RooBifurGauss) RooBifurGauss
int Int_t
Definition: RtypesCore.h:41
RooRealProxy sigmaL
Definition: RooBifurGauss.h:42
STL namespace.
Short_t Abs(Short_t d)
Definition: TMathBase.h:110
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
Double_t Sigma()
Definition: TMath.h:100
RooRealProxy sigmaR
Definition: RooBifurGauss.h:43
RooRealProxy mean
Definition: RooBifurGauss.h:41
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:824
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 atan2(double, double)
#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
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:584
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
Double_t evaluate() const