Logo ROOT   6.08/07
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 \file RooBifurGauss.cxx
20 \class RooBifurGauss
21 \ingroup Roofit
22 
23 Bifurcated Gaussian p.d.f with different widths on left and right
24 side of maximum value
25 **/
26 
27 
28 #include "Riostream.h"
29 #include "TMath.h"
30 #include <math.h>
31 
32 #include "RooBifurGauss.h"
33 #include "RooAbsReal.h"
34 #include "RooMath.h"
35 
36 using namespace std;
37 
39 
40 
41 ////////////////////////////////////////////////////////////////////////////////
42 
43 RooBifurGauss::RooBifurGauss(const char *name, const char *title,
44  RooAbsReal& _x, RooAbsReal& _mean,
45  RooAbsReal& _sigmaL, RooAbsReal& _sigmaR) :
46  RooAbsPdf(name, title),
47  x ("x" , "Dependent" , this, _x),
48  mean ("mean" , "Mean" , this, _mean),
49  sigmaL("sigmaL", "Left Sigma" , this, _sigmaL),
50  sigmaR("sigmaR", "Right Sigma", this, _sigmaR)
51 
52 {
53 }
54 
55 
56 ////////////////////////////////////////////////////////////////////////////////
57 
58 RooBifurGauss::RooBifurGauss(const RooBifurGauss& other, const char* name) :
59  RooAbsPdf(other,name), x("x",this,other.x), mean("mean",this,other.mean),
60  sigmaL("sigmaL",this,other.sigmaL), sigmaR("sigmaR", this, other.sigmaR)
61 {
62 }
63 
64 
65 ////////////////////////////////////////////////////////////////////////////////
66 
68  Double_t arg = x - mean;
69 
70  Double_t coef(0.0);
71 
72  if (arg < 0.0){
73  if (TMath::Abs(sigmaL) > 1e-30) {
74  coef = -0.5/(sigmaL*sigmaL);
75  }
76  } else {
77  if (TMath::Abs(sigmaR) > 1e-30) {
78  coef = -0.5/(sigmaR*sigmaR);
79  }
80  }
81 
82  return exp(coef*arg*arg);
83 }
84 
85 
86 ////////////////////////////////////////////////////////////////////////////////
87 
88 Int_t RooBifurGauss::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
89 {
90  if (matchArgs(allVars,analVars,x)) return 1 ;
91  return 0 ;
92 }
93 
94 
95 ////////////////////////////////////////////////////////////////////////////////
96 
97 Double_t RooBifurGauss::analyticalIntegral(Int_t code, const char* rangeName) const
98 {
99  switch(code) {
100  case 1:
101  {
102  static Double_t root2 = sqrt(2.) ;
103  static Double_t rootPiBy2 = sqrt(atan2(0.0,-1.0)/2.0);
104 
105 // Double_t coefL(0.0), coefR(0.0);
106 // if (TMath::Abs(sigmaL) > 1e-30) {
107 // coefL = -0.5/(sigmaL*sigmaL);
108 // }
109 
110 // if (TMath::Abs(sigmaR) > 1e-30) {
111 // coefR = -0.5/(sigmaR*sigmaR);
112 // }
113 
114  Double_t xscaleL = root2*sigmaL;
115  Double_t xscaleR = root2*sigmaR;
116 
117  Double_t integral = 0.0;
118  if(x.max(rangeName) < mean)
119  {
120  integral = sigmaL * ( RooMath::erf((x.max(rangeName) - mean)/xscaleL) - RooMath::erf((x.min(rangeName) - mean)/xscaleL) );
121  }
122  else if (x.min(rangeName) > mean)
123  {
124  integral = sigmaR * ( RooMath::erf((x.max(rangeName) - mean)/xscaleR) - RooMath::erf((x.min(rangeName) - mean)/xscaleR) );
125  }
126  else
127  {
128  integral = sigmaR*RooMath::erf((x.max(rangeName) - mean)/xscaleR) - sigmaL*RooMath::erf((x.min(rangeName) - mean)/xscaleL);
129  }
130  // return rootPiBy2*(sigmaR*RooMath::erf((x.max(rangeName) - mean)/xscaleR) -
131  // sigmaL*RooMath::erf((x.min(rangeName) - mean)/xscaleL));
132  return integral*rootPiBy2;
133  }
134  }
135 
136  assert(0) ;
137  return 0 ; // to prevent compiler warnings
138 }
RooRealProxy x
Definition: RooBifurGauss.h:40
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
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_t evaluate() const
double sqrt(double)
Double_t x[n]
Definition: legend1.C:17
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
Double_t Sigma()
Definition: TMath.h:100
RooRealProxy sigmaR
Definition: RooBifurGauss.h:43
RooRealProxy mean
Definition: RooBifurGauss.h:41
Bifurcated Gaussian p.d.f with different widths on left and right side of maximum value...
Definition: RooBifurGauss.h:24
Double_t Mean(Long64_t n, const T *a, const Double_t *w=0)
Definition: TMath.h:811
#define ClassImp(name)
Definition: Rtypes.h:279
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 atan2(double, double)
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
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
static std::complex< double > erf(const std::complex< double > z)
complex erf function
Definition: RooMath.cxx:584
double exp(double)
char name[80]
Definition: TGX11.cxx:109