Logo ROOT   6.10/09
Reference Guide
RooGaussModel.h
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * File: $Id: RooGaussModel.h,v 1.21 2007/05/11 09:13:07 verkerke Exp $
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 #ifndef ROO_GAUSS_MODEL
17 #define ROO_GAUSS_MODEL
18 
19 #include "RooResolutionModel.h"
20 #include "RooRealProxy.h"
21 #include "RooMath.h"
22 
23 #include <cmath>
24 #include <complex>
25 
27 public:
28 
38  enum BasisSign { Both=0, Plus=+1, Minus=-1 } ;
39 
40  // Constructors, assignment etc
42  RooGaussModel(const char *name, const char *title, RooRealVar& x,
44  RooGaussModel(const char *name, const char *title, RooRealVar& x,
45  RooAbsReal& mean, RooAbsReal& sigma, RooAbsReal& msSF) ;
46  RooGaussModel(const char *name, const char *title, RooRealVar& x,
47  RooAbsReal& mean, RooAbsReal& sigma, RooAbsReal& meanSF, RooAbsReal& sigmaSF) ;
48  RooGaussModel(const RooGaussModel& other, const char* name=0);
49  virtual TObject* clone(const char* newname) const { return new RooGaussModel(*this,newname) ; }
50  virtual ~RooGaussModel();
51 
52  virtual Int_t basisCode(const char* name) const ;
53  virtual Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
54  virtual Double_t analyticalIntegral(Int_t code, const char* rangeName) const ;
55 
56  Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
57  void generateEvent(Int_t code);
58 
60 
61  void advertiseAymptoticIntegral(Bool_t flag) { _asympInt = flag ; } // added FMV,07/24/03
62 
63 protected:
64 
65  virtual Double_t evaluate() const ;
66  static std::complex<Double_t> evalCerfApprox(Double_t swt, Double_t u, Double_t c);
67 
68  // Calculate exp(-u^2) cwerf(swt*c + i(u+c)), taking care of numerical instabilities
69  static inline std::complex<Double_t> evalCerf(Double_t swt, Double_t u, Double_t c)
70  {
71  std::complex<Double_t> z(swt*c,u+c);
72  return (z.imag()>-4.0) ? (std::exp(-u*u)*RooMath::faddeeva_fast(z)) : evalCerfApprox(swt,u,c);
73  }
74 
75  // Calculate common normalization factors
76  std::complex<Double_t> evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const;
77 
79 
80  Bool_t _asympInt ; // added FMV,07/24/03
81 
86 
87  ClassDef(RooGaussModel,1) // Gaussian Resolution Model
88 };
89 
90 #endif
void advertiseAymptoticIntegral(Bool_t flag)
Definition: RooGaussModel.h:61
RooRealProxy mean
Definition: RooGaussModel.h:82
static std::complex< Double_t > evalCerf(Double_t swt, Double_t u, Double_t c)
Definition: RooGaussModel.h:69
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
static std::complex< Double_t > evalCerfApprox(Double_t swt, Double_t u, Double_t c)
use the approximation: erf(z) = exp(-z*z)/(std::sqrt(pi)*z) to explicitly cancel the divergent exp(y*...
Bool_t _asympInt
Definition: RooGaussModel.h:80
std::complex< Double_t > evalCerfInt(Double_t sign, Double_t wt, Double_t tau, Double_t umin, Double_t umax, Double_t c) const
#define ClassDef(name, id)
Definition: Rtypes.h:297
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t tau
Definition: TRolke.cxx:630
void advertiseFlatScaleFactorIntegral(Bool_t flag)
Definition: RooGaussModel.h:59
virtual Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
virtual Double_t analyticalIntegral(Int_t code, const char *rangeName) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
virtual ~RooGaussModel()
Destructor.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
RooRealProxy ssf
Definition: RooGaussModel.h:85
Bool_t _flatSFInt
Definition: RooGaussModel.h:78
void generateEvent(Int_t code)
Interface for generation of anan event using the algorithm corresponding to the specified code...
static std::complex< double > faddeeva_fast(std::complex< double > z)
evaluate Faddeeva function for complex argument (fast version)
Definition: RooMath.cxx:549
virtual Int_t basisCode(const char *name) const
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Load generatedVars with the subset of directVars that we can generate events for, and return a code t...
const Bool_t kFALSE
Definition: RtypesCore.h:92
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
Mother of all ROOT objects.
Definition: TObject.h:37
you should not use this method at all Int_t Int_t z
Definition: TRolke.cxx:630
virtual Double_t evaluate() const
RooRealProxy msf
Definition: RooGaussModel.h:84
RooRealProxy is the concrete proxy for RooAbsReal objects A RooRealProxy is the general mechanism to ...
Definition: RooRealProxy.h:23
double exp(double)
const Bool_t kTRUE
Definition: RtypesCore.h:91
virtual TObject * clone(const char *newname) const
Definition: RooGaussModel.h:49
RooRealProxy sigma
Definition: RooGaussModel.h:83
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
Definition: RooGaussModel.h:26