00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024 #include "RooFit.h"
00025
00026 #include "Riostream.h"
00027 #include <math.h>
00028
00029 #include "RooMultiVarGaussian.h"
00030 #include "RooAbsReal.h"
00031 #include "RooRealVar.h"
00032 #include "RooRandom.h"
00033 #include "RooMath.h"
00034 #include "RooGlobalFunc.h"
00035 #include "RooConstVar.h"
00036 #include "TDecompChol.h"
00037 #include "RooFitResult.h"
00038
00039 ClassImp(RooMultiVarGaussian)
00040 ;
00041
00042
00043 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
00044 const RooArgList& xvec, const RooArgList& mu, const TMatrixDSym& cov) :
00045 RooAbsPdf(name,title),
00046 _x("x","Observables",this,kTRUE,kFALSE),
00047 _mu("mu","Offset vector",this,kTRUE,kFALSE),
00048 _cov(cov),
00049 _covI(cov),
00050 _z(4)
00051 {
00052 _x.add(xvec) ;
00053
00054 _mu.add(mu) ;
00055
00056 _det = _cov.Determinant() ;
00057
00058
00059 _covI.Invert() ;
00060 }
00061
00062
00063
00064 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
00065 const RooArgList& xvec, const RooFitResult& fr) :
00066 RooAbsPdf(name,title),
00067 _x("x","Observables",this,kTRUE,kFALSE),
00068 _mu("mu","Offset vector",this,kTRUE,kFALSE),
00069 _cov(fr.reducedCovarianceMatrix(xvec)),
00070 _covI(_cov),
00071 _z(4)
00072 {
00073 _det = _cov.Determinant() ;
00074
00075
00076 list<string> munames ;
00077 const RooArgList& fpf = fr.floatParsFinal() ;
00078 for (Int_t i=0 ; i<fpf.getSize() ; i++) {
00079 if (xvec.find(fpf.at(i)->GetName())) {
00080 RooRealVar* parclone = (RooRealVar*) fpf.at(i)->Clone(Form("%s_centralvalue",fpf.at(i)->GetName())) ;
00081 parclone->setConstant(kTRUE) ;
00082 _mu.addOwned(*parclone) ;
00083 munames.push_back(fpf.at(i)->GetName()) ;
00084 }
00085 }
00086
00087
00088 for (list<string>::iterator iter=munames.begin() ; iter!=munames.end() ; iter++) {
00089 RooRealVar* xvar = (RooRealVar*) xvec.find(iter->c_str()) ;
00090 _x.add(*xvar) ;
00091 }
00092
00093
00094 _covI.Invert() ;
00095
00096 }
00097
00098
00099
00100 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
00101 const RooArgList& xvec, const TVectorD& mu, const TMatrixDSym& cov) :
00102 RooAbsPdf(name,title),
00103 _x("x","Observables",this,kTRUE,kFALSE),
00104 _mu("mu","Offset vector",this,kTRUE,kFALSE),
00105 _cov(cov),
00106 _covI(cov),
00107 _z(4)
00108 {
00109 _x.add(xvec) ;
00110
00111 for (Int_t i=0 ; i<mu.GetNrows() ; i++) {
00112 _mu.add(RooFit::RooConst(mu(i))) ;
00113 }
00114
00115 _det = _cov.Determinant() ;
00116
00117
00118 _covI.Invert() ;
00119 }
00120
00121
00122 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
00123 const RooArgList& xvec, const TMatrixDSym& cov) :
00124 RooAbsPdf(name,title),
00125 _x("x","Observables",this,kTRUE,kFALSE),
00126 _mu("mu","Offset vector",this,kTRUE,kFALSE),
00127 _cov(cov),
00128 _covI(cov),
00129 _z(4)
00130 {
00131
00132 _x.add(xvec) ;
00133
00134 for (Int_t i=0 ; i<xvec.getSize() ; i++) {
00135 _mu.add(RooFit::RooConst(0)) ;
00136 }
00137
00138 _det = _cov.Determinant() ;
00139
00140
00141 _covI.Invert() ;
00142 }
00143
00144
00145
00146
00147 RooMultiVarGaussian::RooMultiVarGaussian(const RooMultiVarGaussian& other, const char* name) :
00148 RooAbsPdf(other,name), _x("x",this,other._x), _mu("mu",this,other._mu),
00149 _cov(other._cov), _covI(other._covI), _det(other._det), _z(other._z)
00150 {
00151 }
00152
00153
00154
00155
00156 void RooMultiVarGaussian::syncMuVec() const
00157 {
00158 _muVec.ResizeTo(_mu.getSize()) ;
00159 for (Int_t i=0 ; i<_mu.getSize() ; i++) {
00160 _muVec[i] = ((RooAbsReal*)_mu.at(i))->getVal() ;
00161 }
00162 }
00163
00164
00165
00166 Double_t RooMultiVarGaussian::evaluate() const
00167 {
00168
00169 TVectorD x(_x.getSize()) ;
00170 for (int i=0 ; i<_x.getSize() ; i++) {
00171 x[i] = ((RooAbsReal*)_x.at(i))->getVal() ;
00172 }
00173
00174
00175 syncMuVec() ;
00176 TVectorD x_min_mu = x - _muVec ;
00177
00178 Double_t alpha = x_min_mu * (_covI * x_min_mu) ;
00179 return exp(-0.5*alpha) ;
00180 }
00181
00182
00183
00184
00185 Int_t RooMultiVarGaussian::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
00186 {
00187
00188 if (allVars.getSize()==_x.getSize() && !rangeName) {
00189 analVars.add(allVars) ;
00190 return -1 ;
00191 }
00192
00193 Int_t code(0) ;
00194
00195 Int_t nx = _x.getSize() ;
00196 if (nx>31) {
00197
00198 coutW(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << _x.getSize()
00199 << " observables, analytical integration is only implemented for the first 31 observables" << endl ;
00200 nx=31 ;
00201 }
00202
00203
00204
00205 syncMuVec() ;
00206 for (int i=0 ; i<_x.getSize() ; i++) {
00207
00208 if (allVars.find(_x.at(i)->GetName())) {
00209
00210 RooRealVar* xi = (RooRealVar*)_x.at(i) ;
00211 if (xi->getMin(rangeName)<_muVec(i)-_z*sqrt(_cov(i,i)) && xi->getMax(rangeName) > _muVec(i)+_z*sqrt(_cov(i,i))) {
00212 cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName()
00213 << ") Advertising analytical integral over " << xi->GetName() << " as range is >" << _z << " sigma" << endl ;
00214 code |= (1<<i) ;
00215 analVars.add(*allVars.find(_x.at(i)->GetName())) ;
00216 } else {
00217 cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() << ") Range of " << xi->GetName() << " is <"
00218 << _z << " sigma, relying on numeric integral" << endl ;
00219 }
00220 }
00221 }
00222
00223 return code ;
00224 }
00225
00226
00227
00228
00229 Double_t RooMultiVarGaussian::analyticalIntegral(Int_t code, const char* ) const
00230 {
00231
00232 if (code==-1) {
00233 return pow(2*3.14159268,_x.getSize()/2.)*sqrt(fabs(_det)) ;
00234 }
00235
00236
00237
00238
00239 AnaIntData& aid = anaIntData(code) ;
00240
00241
00242 syncMuVec() ;
00243 TVectorD u(aid.pmap.size()) ;
00244 for (UInt_t i=0 ; i<aid.pmap.size() ; i++) {
00245 u(i) = ((RooAbsReal*)_x.at(aid.pmap[i]))->getVal() - _muVec(aid.pmap[i]) ;
00246 }
00247
00248
00249 Double_t ret = pow(2*3.14159268,aid.nint/2.)/sqrt(fabs(aid.S22det))*exp(-0.5*u*(aid.S22bar*u)) ;
00250
00251 return ret ;
00252 }
00253
00254
00255
00256
00257 RooMultiVarGaussian::AnaIntData& RooMultiVarGaussian::anaIntData(Int_t code) const
00258 {
00259
00260 map<int,AnaIntData>::iterator iter = _anaIntCache.find(code) ;
00261 if (iter != _anaIntCache.end()) {
00262 return iter->second ;
00263 }
00264
00265
00266
00267
00268 vector<int> map1,map2 ;
00269 decodeCode(code,map1,map2) ;
00270
00271
00272
00273
00274 TMatrixDSym S11, S22 ;
00275 TMatrixD S12, S21 ;
00276 blockDecompose(_covI,map1,map2,S11,S12,S21,S22) ;
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291
00292 TMatrixD S22inv(S22) ;
00293 S22inv.Invert() ;
00294 TMatrixD S22bar = S11 - S12*S22inv*S21 ;
00295
00296
00297 AnaIntData& cacheData = _anaIntCache[code] ;
00298 cacheData.S22bar.ResizeTo(S22bar) ;
00299 cacheData.S22bar=S22bar ;
00300 cacheData.S22det= S22.Determinant() ;
00301 cacheData.pmap = map1 ;
00302 cacheData.nint = map2.size() ;
00303
00304 return cacheData ;
00305 }
00306
00307
00308
00309
00310 Int_t RooMultiVarGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t ) const
00311 {
00312
00313 if (directVars.getSize()==_x.getSize()) {
00314 generateVars.add(directVars) ;
00315 return -1 ;
00316 }
00317
00318 Int_t nx = _x.getSize() ;
00319 if (nx>31) {
00320
00321 coutW(Integration) << "RooMultiVarGaussian::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << _x.getSize()
00322 << " observables, partial internal generation is only implemented for the first 31 observables" << endl ;
00323 nx=31 ;
00324 }
00325
00326
00327 Int_t code(0) ;
00328 for (int i=0 ; i<_x.getSize() ; i++) {
00329 RooAbsArg* arg = directVars.find(_x.at(i)->GetName()) ;
00330 if (arg) {
00331 code |= (1<<i) ;
00332 generateVars.add(*arg) ;
00333 }
00334 }
00335
00336 return code ;
00337 }
00338
00339
00340
00341
00342 void RooMultiVarGaussian::initGenerator(Int_t )
00343 {
00344
00345
00346 _genCache.clear() ;
00347
00348 }
00349
00350
00351
00352
00353
00354 void RooMultiVarGaussian::generateEvent(Int_t code)
00355 {
00356
00357 GenData& gd = genData(code) ;
00358 TMatrixD& TU = gd.UT ;
00359 Int_t nobs = TU.GetNcols() ;
00360 vector<int>& omap = gd.omap ;
00361
00362 while(1) {
00363
00364
00365 TVectorD xgen(nobs);
00366 for(Int_t k= 0; k <nobs; k++) {
00367 xgen(k)= RooRandom::gaussian();
00368 }
00369
00370
00371 xgen *= TU ;
00372
00373
00374 if (code == -1) {
00375
00376
00377 xgen += gd.mu1 ;
00378
00379 } else {
00380
00381
00382
00383
00384 TVectorD mubar(gd.mu1) ;
00385 TVectorD x2(gd.pmap.size()) ;
00386 for (UInt_t i=0 ; i<gd.pmap.size() ; i++) {
00387 x2(i) = ((RooAbsReal*)_x.at(gd.pmap[i]))->getVal() ;
00388 }
00389 mubar += gd.S12S22I * (x2 - gd.mu2) ;
00390
00391 xgen += mubar ;
00392
00393 }
00394
00395
00396 Bool_t ok(kTRUE) ;
00397 for (int i=0 ; i<nobs ; i++) {
00398 RooRealVar* xi = (RooRealVar*)_x.at(omap[i]) ;
00399 if (xgen(i)<xi->getMin() || xgen(i)>xi->getMax()) {
00400 ok = kFALSE ;
00401 break ;
00402 } else {
00403 xi->setVal(xgen(i)) ;
00404 }
00405 }
00406
00407
00408
00409 if (ok) {
00410 break ;
00411 }
00412 }
00413
00414 return;
00415 }
00416
00417
00418
00419
00420 RooMultiVarGaussian::GenData& RooMultiVarGaussian::genData(Int_t code) const
00421 {
00422
00423
00424
00425 map<int,GenData>::iterator iter = _genCache.find(code) ;
00426 if (iter != _genCache.end()) {
00427 return iter->second ;
00428 }
00429
00430
00431 GenData& cacheData = _genCache[code] ;
00432
00433 if (code==-1) {
00434
00435
00436 TDecompChol tdc(_cov) ;
00437 tdc.Decompose() ;
00438 TMatrixD U = tdc.GetU() ;
00439 TMatrixD TU(TMatrixD::kTransposed,U) ;
00440
00441
00442 cacheData.UT.ResizeTo(TU) ;
00443 cacheData.UT = TU ;
00444 cacheData.omap.resize(_x.getSize()) ;
00445 for (int i=0 ; i<_x.getSize() ; i++) {
00446 cacheData.omap[i] = i ;
00447 }
00448 syncMuVec() ;
00449 cacheData.mu1.ResizeTo(_muVec) ;
00450 cacheData.mu1 = _muVec ;
00451
00452 } else {
00453
00454
00455 vector<int> map1, map2 ;
00456 decodeCode(code,map2,map1) ;
00457
00458
00459 TMatrixDSym S11, S22 ;
00460 TMatrixD S12, S21 ;
00461 blockDecompose(_cov,map1,map2,S11,S12,S21,S22) ;
00462
00463
00464
00465
00466
00467
00468
00469
00470 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
00471 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
00472
00473
00474 TDecompChol tdc(S22bar) ;
00475 tdc.Decompose() ;
00476 TMatrixD U = tdc.GetU() ;
00477 TMatrixD TU(TMatrixD::kTransposed,U) ;
00478
00479
00480 TVectorD mu1(map1.size()),mu2(map2.size()) ;
00481 syncMuVec() ;
00482 for (UInt_t i=0 ; i<map1.size() ; i++) {
00483 mu1(i) = _muVec(map1[i]) ;
00484 }
00485 for (UInt_t i=0 ; i<map2.size() ; i++) {
00486 mu2(i) = _muVec(map2[i]) ;
00487 }
00488
00489
00490 TMatrixD S12S22Inv = S12 * S22Inv ;
00491
00492
00493 cacheData.UT.ResizeTo(TU) ;
00494 cacheData.UT = TU ;
00495 cacheData.omap = map1 ;
00496 cacheData.pmap = map2 ;
00497 cacheData.mu1.ResizeTo(mu1) ;
00498 cacheData.mu2.ResizeTo(mu2) ;
00499 cacheData.mu1 = mu1 ;
00500 cacheData.mu2 = mu2 ;
00501 cacheData.S12S22I.ResizeTo(S12S22Inv) ;
00502 cacheData.S12S22I = S12S22Inv ;
00503
00504 }
00505
00506
00507 return cacheData ;
00508 }
00509
00510
00511
00512
00513
00514 void RooMultiVarGaussian::decodeCode(Int_t code, vector<int>& map1, vector<int>& map2) const
00515 {
00516
00517
00518 map1.clear() ;
00519 map2.clear() ;
00520 for (int i=0 ; i<_x.getSize() ; i++) {
00521 if (code & (1<<i)) {
00522 map2.push_back(i) ;
00523 } else {
00524 map1.push_back(i) ;
00525 }
00526 }
00527 }
00528
00529
00530
00531 void RooMultiVarGaussian::blockDecompose(const TMatrixD& input, const vector<int>& map1, const vector<int>& map2, TMatrixDSym& S11, TMatrixD& S12, TMatrixD& S21, TMatrixDSym& S22)
00532 {
00533
00534
00535
00536
00537 S11.ResizeTo(map1.size(),map1.size()) ;
00538 S12.ResizeTo(map1.size(),map2.size()) ;
00539 S21.ResizeTo(map2.size(),map1.size()) ;
00540 S22.ResizeTo(map2.size(),map2.size()) ;
00541
00542 for (UInt_t i=0 ; i<map1.size() ; i++) {
00543 for (UInt_t j=0 ; j<map1.size() ; j++)
00544 S11(i,j) = input(map1[i],map1[j]) ;
00545 for (UInt_t j=0 ; j<map2.size() ; j++)
00546 S12(i,j) = input(map1[i],map2[j]) ;
00547 }
00548 for (UInt_t i=0 ; i<map2.size() ; i++) {
00549 for (UInt_t j=0 ; j<map1.size() ; j++)
00550 S21(i,j) = input(map2[i],map1[j]) ;
00551 for (UInt_t j=0 ; j<map2.size() ; j++)
00552 S22(i,j) = input(map2[i],map2[j]) ;
00553 }
00554
00555 }
00556
00557
00558
00559
00560
00561
00562