Logo ROOT   6.10/09
Reference Guide
RooMultiVarGaussian.h
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * File: $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 #ifndef ROO_MULTI_VAR_GAUSSIAN
17 #define ROO_MULTI_VAR_GAUSSIAN
18 
19 #include "RooAbsPdf.h"
20 #include "RooListProxy.h"
21 #include "TMatrixDSym.h"
22 #include "TMatrixD.h"
23 #include "TVectorD.h"
24 
25 class RooRealVar;
26 class RooFitResult ;
27 
28 #include <map>
29 #include <vector>
30 
32 public:
33 
35  RooMultiVarGaussian(const char *name, const char *title, const RooArgList& xvec, const RooArgList& mu, const TMatrixDSym& covMatrix) ;
36  RooMultiVarGaussian(const char *name, const char *title, const RooArgList& xvec, const RooFitResult& fr, Bool_t reduceToConditional=kTRUE) ;
37  RooMultiVarGaussian(const char *name, const char *title, const RooArgList& xvec, const TVectorD& mu, const TMatrixDSym& covMatrix) ;
38  RooMultiVarGaussian(const char *name, const char *title, const RooArgList& xvec,const TMatrixDSym& covMatrix) ;
39  void setAnaIntZ(Double_t z) { _z = z ; }
40 
41  RooMultiVarGaussian(const RooMultiVarGaussian& other, const char* name=0) ;
42  virtual TObject* clone(const char* newname) const { return new RooMultiVarGaussian(*this,newname); }
43  inline virtual ~RooMultiVarGaussian() { }
44 
45  Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const ;
46  Double_t analyticalIntegral(Int_t code, const char* rangeName=0) const ;
47 
48  Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const;
49  void initGenerator(Int_t code) ;
50  void generateEvent(Int_t code);
51 
52  const TMatrixDSym& covarianceMatrix() const { return _cov ; }
53 
54  class AnaIntData {
55  public:
58  std::vector<int> pmap ;
60  } ;
61 
62  class GenData {
63  public:
65  std::vector<int> omap ;
66  std::vector<int> pmap ;
70  } ;
71 
72  class BitBlock {
73  public:
74  BitBlock() : b0(0), b1(0), b2(0), b3(0) {} ;
75 
76  void setBit(Int_t ibit) ;
77  Bool_t getBit(Int_t ibit) ;
78  Bool_t operator==(const BitBlock& other) ;
79 
84  } ;
85 
86  static void blockDecompose(const TMatrixD& input, const std::vector<int>& map1, const std::vector<int>& map2, TMatrixDSym& S11, TMatrixD& S12, TMatrixD& S21, TMatrixDSym& S22) ;
87 
88 protected:
89 
90  void decodeCode(Int_t code, std::vector<int>& map1, std::vector<int>& map2) const;
91  AnaIntData& anaIntData(Int_t code) const ;
92  GenData& genData(Int_t code) const ;
93 
94  mutable std::map<int,AnaIntData> _anaIntCache ; //!
95  mutable std::map<int,GenData> _genCache ; //!
96 
97  mutable std::vector<BitBlock> _aicMap ; //!
98 
105 
106  void syncMuVec() const ;
107  mutable TVectorD _muVec ; //! Do not persist
108 
109  Double_t evaluate() const ;
110 
111 private:
112 
113  ClassDef(RooMultiVarGaussian,1) // Multivariate Gaussian PDF with correlations
114 };
115 
116 #endif
std::map< int, GenData > _genCache
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Special case: generate all observables.
void setAnaIntZ(Double_t z)
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
static void blockDecompose(const TMatrixD &input, const std::vector< int > &map1, const std::vector< int > &map2, TMatrixDSym &S11, TMatrixD &S12, TMatrixD &S21, TMatrixDSym &S22)
Block decomposition of covI according to given maps of observables.
const TMatrixDSym & covarianceMatrix() const
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
void decodeCode(Int_t code, std::vector< int > &map1, std::vector< int > &map2) const
Decode analytical integration/generation code into index map of integrated/generated (map2) and non-i...
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Handle full integral here.
Double_t evaluate() const
Do not persist.
std::vector< BitBlock > _aicMap
std::map< int, AnaIntData > _anaIntCache
#define ClassDef(name, id)
Definition: Rtypes.h:297
void initGenerator(Int_t code)
Clear the GenData cache as its content is not invariant under changes in the mu vector.
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void generateEvent(Int_t code)
Retrieve generator config from cache.
GenData & genData(Int_t code) const
WVE – CHECK THAT GENDATA IS VALID GIVEN CURRENT VALUES OF _MU.
RooListProxy is the concrete proxy for RooArgList objects.
Definition: RooListProxy.h:25
Multivariate Gaussian p.d.f.
Bool_t operator==(Double_t value) const
Equality operator comparing to a Double_t.
Definition: RooAbsReal.cxx:194
double Double_t
Definition: RtypesCore.h:55
AnaIntData & anaIntData(Int_t code) const
Check if cache entry was previously created.
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
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual TObject * clone(const char *newname) const
const Bool_t kTRUE
Definition: RtypesCore.h:91