42 ClassImp(RooMultiVarGaussian);
47 RooMultiVarGaussian::RooMultiVarGaussian(
const char *name,
const char *title,
48 const RooArgList& xvec,
const RooArgList& mu,
const TMatrixDSym& cov) :
49 RooAbsPdf(name,title),
50 _x(
"x",
"Observables",this,kTRUE,kFALSE),
51 _mu(
"mu",
"Offset vector",this,kTRUE,kFALSE),
60 _det = _cov.Determinant() ;
69 RooMultiVarGaussian::RooMultiVarGaussian(
const char *name,
const char *title,
70 const RooArgList& xvec,
const RooFitResult& fr, Bool_t reduceToConditional) :
71 RooAbsPdf(name,title),
72 _x(
"x",
"Observables",this,kTRUE,kFALSE),
73 _mu(
"mu",
"Offset vector",this,kTRUE,kFALSE),
74 _cov(reduceToConditional ? fr.conditionalCovarianceMatrix(xvec) : fr.reducedCovarianceMatrix(xvec)),
78 _det = _cov.Determinant() ;
81 list<string> munames ;
82 const RooArgList& fpf = fr.floatParsFinal() ;
83 for (Int_t i=0 ; i<fpf.getSize() ; i++) {
84 if (xvec.find(fpf.at(i)->GetName())) {
85 RooRealVar* parclone = (RooRealVar*) fpf.at(i)->Clone(Form(
"%s_centralvalue",fpf.at(i)->GetName())) ;
86 parclone->setConstant(kTRUE) ;
87 _mu.addOwned(*parclone) ;
88 munames.push_back(fpf.at(i)->GetName()) ;
93 for (list<string>::iterator iter=munames.begin() ; iter!=munames.end() ; ++iter) {
94 RooRealVar* xvar = (RooRealVar*) xvec.find(iter->c_str()) ;
106 RooMultiVarGaussian::RooMultiVarGaussian(
const char *name,
const char *title,
107 const RooArgList& xvec,
const TVectorD& mu,
const TMatrixDSym& cov) :
108 RooAbsPdf(name,title),
109 _x(
"x",
"Observables",this,kTRUE,kFALSE),
110 _mu(
"mu",
"Offset vector",this,kTRUE,kFALSE),
117 for (Int_t i=0 ; i<mu.GetNrows() ; i++) {
118 _mu.add(RooFit::RooConst(mu(i))) ;
121 _det = _cov.Determinant() ;
129 RooMultiVarGaussian::RooMultiVarGaussian(
const char *name,
const char *title,
130 const RooArgList& xvec,
const TMatrixDSym& cov) :
131 RooAbsPdf(name,title),
132 _x(
"x",
"Observables",this,kTRUE,kFALSE),
133 _mu(
"mu",
"Offset vector",this,kTRUE,kFALSE),
140 for (Int_t i=0 ; i<xvec.getSize() ; i++) {
141 _mu.add(RooFit::RooConst(0)) ;
144 _det = _cov.Determinant() ;
154 RooMultiVarGaussian::RooMultiVarGaussian(
const RooMultiVarGaussian& other,
const char* name) :
155 RooAbsPdf(other,name), _aicMap(other._aicMap), _x(
"x",this,other._x), _mu(
"mu",this,other._mu),
156 _cov(other._cov), _covI(other._covI), _det(other._det), _z(other._z)
164 void RooMultiVarGaussian::syncMuVec()
const
166 _muVec.ResizeTo(_mu.getSize()) ;
167 for (Int_t i=0 ; i<_mu.getSize() ; i++) {
168 _muVec[i] = ((RooAbsReal*)_mu.at(i))->getVal() ;
176 Double_t RooMultiVarGaussian::evaluate()
const
178 TVectorD x(_x.getSize()) ;
179 for (
int i=0 ; i<_x.getSize() ; i++) {
180 x[i] = ((RooAbsReal*)_x.at(i))->getVal() ;
185 TVectorD x_min_mu = x - _muVec ;
187 Double_t alpha = x_min_mu * (_covI * x_min_mu) ;
188 return exp(-0.5*alpha) ;
195 Int_t RooMultiVarGaussian::getAnalyticalIntegral(RooArgSet& allVarsIn, RooArgSet& analVars,
const char* rangeName)
const
197 RooArgSet allVars(allVarsIn) ;
200 for (Int_t i=0 ; i<_x.getSize() ; i++) {
201 if (allVars.contains(*_x.at(i))) {
202 allVars.remove(*_mu.at(i),kTRUE,kTRUE) ;
208 if (allVars.getSize()==_x.getSize() && !rangeName) {
209 analVars.add(allVars) ;
215 Int_t nx = _x.getSize() ;
218 coutW(Integration) <<
"RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() <<
") WARNING: p.d.f. has " << _x.getSize()
219 <<
" observables, analytical integration is only implemented for the first 127 observables" << endl ;
226 Bool_t anyBits(kFALSE) ;
228 for (
int i=0 ; i<_x.getSize() ; i++) {
231 if (allVars.find(_x.at(i)->GetName())) {
233 RooRealVar* xi = (RooRealVar*)_x.at(i) ;
234 if (xi->getMin(rangeName)<_muVec(i)-_z*sqrt(_cov(i,i)) && xi->getMax(rangeName) > _muVec(i)+_z*sqrt(_cov(i,i))) {
235 cxcoutD(Integration) <<
"RooMultiVarGaussian::getAnalyticalIntegral(" << GetName()
236 <<
") Advertising analytical integral over " << xi->GetName() <<
" as range is >" << _z <<
" sigma" << endl ;
239 analVars.add(*allVars.find(_x.at(i)->GetName())) ;
241 cxcoutD(Integration) <<
"RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() <<
") Range of " << xi->GetName() <<
" is <"
242 << _z <<
" sigma, relying on numeric integral" << endl ;
247 if (allVars.find(_mu.at(i)->GetName())) {
249 RooRealVar* pi = (RooRealVar*)_mu.at(i) ;
250 if (pi->getMin(rangeName)<_muVec(i)-_z*sqrt(_cov(i,i)) && pi->getMax(rangeName) > _muVec(i)+_z*sqrt(_cov(i,i))) {
251 cxcoutD(Integration) <<
"RooMultiVarGaussian::getAnalyticalIntegral(" << GetName()
252 <<
") Advertising analytical integral over " << pi->GetName() <<
" as range is >" << _z <<
" sigma" << endl ;
255 analVars.add(*allVars.find(_mu.at(i)->GetName())) ;
257 cxcoutD(Integration) <<
"RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() <<
") Range of " << pi->GetName() <<
" is <"
258 << _z <<
" sigma, relying on numeric integral" << endl ;
271 for (UInt_t i=0 ; i<_aicMap.size() ; i++) {
272 if (_aicMap[i]==bits) {
277 _aicMap.push_back(bits) ;
278 code = _aicMap.size() ;
289 Double_t RooMultiVarGaussian::analyticalIntegral(Int_t code,
const char* )
const
292 return pow(2*3.14159268,_x.getSize()/2.)*sqrt(fabs(_det)) ;
298 AnaIntData& aid = anaIntData(code) ;
302 TVectorD u(aid.pmap.size()) ;
303 for (UInt_t i=0 ; i<aid.pmap.size() ; i++) {
304 u(i) = ((RooAbsReal*)_x.at(aid.pmap[i]))->getVal() - _muVec(aid.pmap[i]) ;
308 Double_t ret = pow(2*3.14159268,aid.nint/2.)/sqrt(fabs(aid.S22det))*exp(-0.5*u*(aid.S22bar*u)) ;
318 RooMultiVarGaussian::AnaIntData& RooMultiVarGaussian::anaIntData(Int_t code)
const
320 map<int,AnaIntData>::iterator iter = _anaIntCache.find(code) ;
321 if (iter != _anaIntCache.end()) {
322 return iter->second ;
328 vector<int> map1,map2 ;
329 decodeCode(code,map1,map2) ;
334 TMatrixDSym S11, S22 ;
336 blockDecompose(_covI,map1,map2,S11,S12,S21,S22) ;
352 TMatrixD S22inv(S22) ;
354 TMatrixD S22bar = S11 - S12*S22inv*S21 ;
357 AnaIntData& cacheData = _anaIntCache[code] ;
358 cacheData.S22bar.ResizeTo(S22bar) ;
359 cacheData.S22bar=S22bar ;
360 cacheData.S22det= S22.Determinant() ;
361 cacheData.pmap = map1 ;
362 cacheData.nint = map2.size() ;
372 Int_t RooMultiVarGaussian::getGenerator(
const RooArgSet& directVars, RooArgSet &generateVars, Bool_t )
const
374 if (directVars.getSize()==_x.getSize()) {
375 generateVars.add(directVars) ;
379 Int_t nx = _x.getSize() ;
382 coutW(Integration) <<
"RooMultiVarGaussian::getGenerator(" << GetName() <<
") WARNING: p.d.f. has " << _x.getSize()
383 <<
" observables, partial internal generation is only implemented for the first 127 observables" << endl ;
390 for (
int i=0 ; i<_x.getSize() ; i++) {
391 RooAbsArg* arg = directVars.find(_x.at(i)->GetName()) ;
395 generateVars.add(*arg) ;
400 for (UInt_t i=0 ; i<_aicMap.size() ; i++) {
401 if (_aicMap[i]==bits) {
406 _aicMap.push_back(bits) ;
407 code = _aicMap.size() ;
420 void RooMultiVarGaussian::initGenerator(Int_t )
432 void RooMultiVarGaussian::generateEvent(Int_t code)
434 GenData& gd = genData(code) ;
435 TMatrixD& TU = gd.UT ;
436 Int_t nobs = TU.GetNcols() ;
437 vector<int>& omap = gd.omap ;
443 for(Int_t k= 0; k <nobs; k++) {
444 xgen(k)= RooRandom::gaussian();
461 TVectorD mubar(gd.mu1) ;
462 TVectorD x2(gd.pmap.size()) ;
463 for (UInt_t i=0 ; i<gd.pmap.size() ; i++) {
464 x2(i) = ((RooAbsReal*)_x.at(gd.pmap[i]))->getVal() ;
466 mubar += gd.S12S22I * (x2 - gd.mu2) ;
474 for (
int i=0 ; i<nobs ; i++) {
475 RooRealVar* xi = (RooRealVar*)_x.at(omap[i]) ;
476 if (xgen(i)<xi->getMin() || xgen(i)>xi->getMax()) {
480 xi->setVal(xgen(i)) ;
499 RooMultiVarGaussian::GenData& RooMultiVarGaussian::genData(Int_t code)
const
502 map<int,GenData>::iterator iter = _genCache.find(code) ;
503 if (iter != _genCache.end()) {
504 return iter->second ;
508 GenData& cacheData = _genCache[code] ;
513 TDecompChol tdc(_cov) ;
515 TMatrixD U = tdc.GetU() ;
516 TMatrixD TU(TMatrixD::kTransposed,U) ;
519 cacheData.UT.ResizeTo(TU) ;
521 cacheData.omap.resize(_x.getSize()) ;
522 for (
int i=0 ; i<_x.getSize() ; i++) {
523 cacheData.omap[i] = i ;
526 cacheData.mu1.ResizeTo(_muVec) ;
527 cacheData.mu1 = _muVec ;
532 vector<int> map1, map2 ;
533 decodeCode(code,map2,map1) ;
536 TMatrixDSym S11, S22 ;
538 blockDecompose(_cov,map1,map2,S11,S12,S21,S22) ;
547 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
548 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
551 TDecompChol tdc(S22bar) ;
553 TMatrixD U = tdc.GetU() ;
554 TMatrixD TU(TMatrixD::kTransposed,U) ;
557 TVectorD mu1(map1.size()),mu2(map2.size()) ;
559 for (UInt_t i=0 ; i<map1.size() ; i++) {
560 mu1(i) = _muVec(map1[i]) ;
562 for (UInt_t i=0 ; i<map2.size() ; i++) {
563 mu2(i) = _muVec(map2[i]) ;
567 TMatrixD S12S22Inv = S12 * S22Inv ;
570 cacheData.UT.ResizeTo(TU) ;
572 cacheData.omap = map1 ;
573 cacheData.pmap = map2 ;
574 cacheData.mu1.ResizeTo(mu1) ;
575 cacheData.mu2.ResizeTo(mu2) ;
576 cacheData.mu1 = mu1 ;
577 cacheData.mu2 = mu2 ;
578 cacheData.S12S22I.ResizeTo(S12S22Inv) ;
579 cacheData.S12S22I = S12S22Inv ;
594 void RooMultiVarGaussian::decodeCode(Int_t code, vector<int>& map1, vector<int>& map2)
const
596 if (code<0 || code> (Int_t)_aicMap.size()) {
597 cout <<
"RooMultiVarGaussian::decodeCode(" << GetName() <<
") ERROR don't have bit pattern for code " << code << endl ;
598 throw string(
"RooMultiVarGaussian::decodeCode() ERROR don't have bit pattern for code") ;
601 BitBlock b = _aicMap[code-1] ;
604 for (
int i=0 ; i<_x.getSize() ; i++) {
617 void RooMultiVarGaussian::blockDecompose(
const TMatrixD& input,
const vector<int>& map1,
const vector<int>& map2, TMatrixDSym& S11, TMatrixD& S12, TMatrixD& S21, TMatrixDSym& S22)
621 S11.ResizeTo(map1.size(),map1.size()) ;
622 S12.ResizeTo(map1.size(),map2.size()) ;
623 S21.ResizeTo(map2.size(),map1.size()) ;
624 S22.ResizeTo(map2.size(),map2.size()) ;
626 for (UInt_t i=0 ; i<map1.size() ; i++) {
627 for (UInt_t j=0 ; j<map1.size() ; j++)
628 S11(i,j) = input(map1[i],map1[j]) ;
629 for (UInt_t j=0 ; j<map2.size() ; j++)
630 S12(i,j) = input(map1[i],map2[j]) ;
632 for (UInt_t i=0 ; i<map2.size() ; i++) {
633 for (UInt_t j=0 ; j<map1.size() ; j++)
634 S21(i,j) = input(map2[i],map1[j]) ;
635 for (UInt_t j=0 ; j<map2.size() ; j++)
636 S22(i,j) = input(map2[i],map2[j]) ;
642 void RooMultiVarGaussian::BitBlock::setBit(Int_t ibit)
644 if (ibit<32) { b0 |= (1<<ibit) ; return ; }
645 if (ibit<64) { b1 |= (1<<(ibit-32)) ; return ; }
646 if (ibit<96) { b2 |= (1<<(ibit-64)) ; return ; }
647 if (ibit<128) { b3 |= (1<<(ibit-96)) ; return ; }
650 Bool_t RooMultiVarGaussian::BitBlock::getBit(Int_t ibit)
652 if (ibit<32)
return (b0 & (1<<ibit)) ;
653 if (ibit<64)
return (b1 & (1<<(ibit-32))) ;
654 if (ibit<96)
return (b2 & (1<<(ibit-64))) ;
655 if (ibit<128)
return (b3 & (1<<(ibit-96))) ;
659 Bool_t RooMultiVarGaussian::BitBlock::operator==(
const BitBlock& other)
661 if (b0 != other.b0)
return kFALSE ;
662 if (b1 != other.b1)
return kFALSE ;
663 if (b2 != other.b2)
return kFALSE ;
664 if (b3 != other.b3)
return kFALSE ;