ROOT  6.06/09
Reference Guide
RooMultiVarGaussian.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$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 
17 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // Multivariate Gaussian p.d.f. with correlations
21 // END_HTML
22 //
23 
24 #include "RooFit.h"
25 
26 #include "Riostream.h"
27 #include <math.h>
28 
29 #include "RooMultiVarGaussian.h"
30 #include "RooAbsReal.h"
31 #include "RooRealVar.h"
32 #include "RooRandom.h"
33 #include "RooMath.h"
34 #include "RooGlobalFunc.h"
35 #include "RooConstVar.h"
36 #include "TDecompChol.h"
37 #include "RooFitResult.h"
38 
39 using namespace std;
40 
42  ;
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 
46 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
47  const RooArgList& xvec, const RooArgList& mu, const TMatrixDSym& cov) :
48  RooAbsPdf(name,title),
49  _x("x","Observables",this,kTRUE,kFALSE),
50  _mu("mu","Offset vector",this,kTRUE,kFALSE),
51  _cov(cov),
52  _covI(cov),
53  _z(4)
54 {
55  _x.add(xvec) ;
56 
57  _mu.add(mu) ;
58 
59  _det = _cov.Determinant() ;
60 
61  // Invert covariance matrix
62  _covI.Invert() ;
63 }
64 
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 
68 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
69  const RooArgList& xvec, const RooFitResult& fr, Bool_t reduceToConditional) :
70  RooAbsPdf(name,title),
71  _x("x","Observables",this,kTRUE,kFALSE),
72  _mu("mu","Offset vector",this,kTRUE,kFALSE),
73  _cov(reduceToConditional ? fr.conditionalCovarianceMatrix(xvec) : fr.reducedCovarianceMatrix(xvec)),
74  _covI(_cov),
75  _z(4)
76 {
77  _det = _cov.Determinant() ;
78 
79  // Fill mu vector with constant RooRealVars
80  list<string> munames ;
81  const RooArgList& fpf = fr.floatParsFinal() ;
82  for (Int_t i=0 ; i<fpf.getSize() ; i++) {
83  if (xvec.find(fpf.at(i)->GetName())) {
84  RooRealVar* parclone = (RooRealVar*) fpf.at(i)->Clone(Form("%s_centralvalue",fpf.at(i)->GetName())) ;
85  parclone->setConstant(kTRUE) ;
86  _mu.addOwned(*parclone) ;
87  munames.push_back(fpf.at(i)->GetName()) ;
88  }
89  }
90 
91  // Fill X vector in same order as mu vector
92  for (list<string>::iterator iter=munames.begin() ; iter!=munames.end() ; iter++) {
93  RooRealVar* xvar = (RooRealVar*) xvec.find(iter->c_str()) ;
94  _x.add(*xvar) ;
95  }
96 
97  // Invert covariance matrix
98  _covI.Invert() ;
99 
100 }
101 
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 
105 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
106  const RooArgList& xvec, const TVectorD& mu, const TMatrixDSym& cov) :
107  RooAbsPdf(name,title),
108  _x("x","Observables",this,kTRUE,kFALSE),
109  _mu("mu","Offset vector",this,kTRUE,kFALSE),
110  _cov(cov),
111  _covI(cov),
112  _z(4)
113 {
114  _x.add(xvec) ;
115 
116  for (Int_t i=0 ; i<mu.GetNrows() ; i++) {
117  _mu.add(RooFit::RooConst(mu(i))) ;
118  }
119 
120  _det = _cov.Determinant() ;
121 
122  // Invert covariance matrix
123  _covI.Invert() ;
124 }
125 
126 ////////////////////////////////////////////////////////////////////////////////
127 
128 RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
129  const RooArgList& xvec, const TMatrixDSym& cov) :
130  RooAbsPdf(name,title),
131  _x("x","Observables",this,kTRUE,kFALSE),
132  _mu("mu","Offset vector",this,kTRUE,kFALSE),
133  _cov(cov),
134  _covI(cov),
135  _z(4)
136 {
137  _x.add(xvec) ;
138 
139  for (Int_t i=0 ; i<xvec.getSize() ; i++) {
140  _mu.add(RooFit::RooConst(0)) ;
141  }
142 
143  _det = _cov.Determinant() ;
144 
145  // Invert covariance matrix
146  _covI.Invert() ;
147 }
148 
149 
150 
151 ////////////////////////////////////////////////////////////////////////////////
152 
154  RooAbsPdf(other,name), _aicMap(other._aicMap), _x("x",this,other._x), _mu("mu",this,other._mu),
155  _cov(other._cov), _covI(other._covI), _det(other._det), _z(other._z)
156 {
157 }
158 
159 
160 
161 ////////////////////////////////////////////////////////////////////////////////
162 
164 {
166  for (Int_t i=0 ; i<_mu.getSize() ; i++) {
167  _muVec[i] = ((RooAbsReal*)_mu.at(i))->getVal() ;
168  }
169 }
170 
171 
172 ////////////////////////////////////////////////////////////////////////////////
173 /// Represent observables as vector
174 
176 {
177  TVectorD x(_x.getSize()) ;
178  for (int i=0 ; i<_x.getSize() ; i++) {
179  x[i] = ((RooAbsReal*)_x.at(i))->getVal() ;
180  }
181 
182  // Calculate return value
183  syncMuVec() ;
184  TVectorD x_min_mu = x - _muVec ;
185 
186  Double_t alpha = x_min_mu * (_covI * x_min_mu) ;
187  return exp(-0.5*alpha) ;
188 }
189 
190 
191 
192 ////////////////////////////////////////////////////////////////////////////////
193 
194 Int_t RooMultiVarGaussian::getAnalyticalIntegral(RooArgSet& allVarsIn, RooArgSet& analVars, const char* rangeName) const
195 {
196  RooArgSet allVars(allVarsIn) ;
197 
198  // If allVars contains x_i it cannot contain mu_i
199  for (Int_t i=0 ; i<_x.getSize() ; i++) {
200  if (allVars.contains(*_x.at(i))) {
201  allVars.remove(*_mu.at(i),kTRUE,kTRUE) ;
202  }
203  }
204 
205 
206  // Analytical integral known over all observables
207  if (allVars.getSize()==_x.getSize() && !rangeName) {
208  analVars.add(allVars) ;
209  return -1 ;
210  }
211 
212  Int_t code(0) ;
213 
214  Int_t nx = _x.getSize() ;
215  if (nx>127) {
216  // Warn that analytical integration is only provided for the first 127 observables
217  coutW(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() << ") WARNING: p.d.f. has " << _x.getSize()
218  << " observables, analytical integration is only implemented for the first 127 observables" << endl ;
219  nx=127 ;
220  }
221 
222  // Advertise partial analytical integral over all observables for which is wide enough to
223  // use asymptotic integral calculation
224  BitBlock bits ;
225  Bool_t anyBits(kFALSE) ;
226  syncMuVec() ;
227  for (int i=0 ; i<_x.getSize() ; i++) {
228 
229  // Check if integration over observable #i is requested
230  if (allVars.find(_x.at(i)->GetName())) {
231  // Check if range is wider than Z sigma
232  RooRealVar* xi = (RooRealVar*)_x.at(i) ;
233  if (xi->getMin(rangeName)<_muVec(i)-_z*sqrt(_cov(i,i)) && xi->getMax(rangeName) > _muVec(i)+_z*sqrt(_cov(i,i))) {
234  cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName()
235  << ") Advertising analytical integral over " << xi->GetName() << " as range is >" << _z << " sigma" << endl ;
236  bits.setBit(i) ;
237  anyBits = kTRUE ;
238  analVars.add(*allVars.find(_x.at(i)->GetName())) ;
239  } else {
240  cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() << ") Range of " << xi->GetName() << " is <"
241  << _z << " sigma, relying on numeric integral" << endl ;
242  }
243  }
244 
245  // Check if integration over parameter #i is requested
246  if (allVars.find(_mu.at(i)->GetName())) {
247  // Check if range is wider than Z sigma
248  RooRealVar* pi = (RooRealVar*)_mu.at(i) ;
249  if (pi->getMin(rangeName)<_muVec(i)-_z*sqrt(_cov(i,i)) && pi->getMax(rangeName) > _muVec(i)+_z*sqrt(_cov(i,i))) {
250  cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName()
251  << ") Advertising analytical integral over " << pi->GetName() << " as range is >" << _z << " sigma" << endl ;
252  bits.setBit(i) ;
253  anyBits = kTRUE ;
254  analVars.add(*allVars.find(_mu.at(i)->GetName())) ;
255  } else {
256  cxcoutD(Integration) << "RooMultiVarGaussian::getAnalyticalIntegral(" << GetName() << ") Range of " << pi->GetName() << " is <"
257  << _z << " sigma, relying on numeric integral" << endl ;
258  }
259  }
260 
261 
262  }
263 
264  // Full numeric integration over requested observables maps always to code zero
265  if (!anyBits) {
266  return 0 ;
267  }
268 
269  // Map BitBlock into return code
270  for (UInt_t i=0 ; i<_aicMap.size() ; i++) {
271  if (_aicMap[i]==bits) {
272  code = i+1 ;
273  }
274  }
275  if (code==0) {
276  _aicMap.push_back(bits) ;
277  code = _aicMap.size() ;
278  }
279 
280  return code ;
281 }
282 
283 
284 
285 ////////////////////////////////////////////////////////////////////////////////
286 /// Handle full integral here
287 
288 Double_t RooMultiVarGaussian::analyticalIntegral(Int_t code, const char* /*rangeName*/) const
289 {
290  if (code==-1) {
291  return pow(2*3.14159268,_x.getSize()/2.)*sqrt(fabs(_det)) ;
292  }
293 
294  // Handle partial integrals here
295 
296  // Retrieve |S22|, S22bar from cache
297  AnaIntData& aid = anaIntData(code) ;
298 
299  // Fill position vector for non-integrated observables
300  syncMuVec() ;
301  TVectorD u(aid.pmap.size()) ;
302  for (UInt_t i=0 ; i<aid.pmap.size() ; i++) {
303  u(i) = ((RooAbsReal*)_x.at(aid.pmap[i]))->getVal() - _muVec(aid.pmap[i]) ;
304  }
305 
306  // Calculate partial integral
307  Double_t ret = pow(2*3.14159268,aid.nint/2.)/sqrt(fabs(aid.S22det))*exp(-0.5*u*(aid.S22bar*u)) ;
308 
309  return ret ;
310 }
311 
312 
313 
314 ////////////////////////////////////////////////////////////////////////////////
315 /// Check if cache entry was previously created
316 
318 {
319  map<int,AnaIntData>::iterator iter = _anaIntCache.find(code) ;
320  if (iter != _anaIntCache.end()) {
321  return iter->second ;
322  }
323 
324  // Calculate cache contents
325 
326  // Decode integration code
327  vector<int> map1,map2 ;
328  decodeCode(code,map1,map2) ;
329 
330  // Rearrage observables so that all non-integrated observables
331  // go first (preserving relative order) and all integrated observables
332  // go last (preserving relative order)
333  TMatrixDSym S11, S22 ;
334  TMatrixD S12, S21 ;
335  blockDecompose(_covI,map1,map2,S11,S12,S21,S22) ;
336 
337  // Begin calculation of partial integrals
338  // ___
339  // sqrt(2pi)^(#intObs) (-0.5 * u1T S22 u1 )
340  // I = ------------------- * e
341  // sqrt(|det(S22)|)
342  // ___
343  // Where S22 is the sub-matrix of covI for the integrated observables and S22
344  // is the Schur complement of S22
345  // ___ -1
346  // S22 = S11 - S12 * S22 * S21
347  //
348  // and u1 is the vector of non-integrated observables
349 
350  // Calculate Schur complement S22bar
351  TMatrixD S22inv(S22) ;
352  S22inv.Invert() ;
353  TMatrixD S22bar = S11 - S12*S22inv*S21 ;
354 
355  // Create new cache entry
356  AnaIntData& cacheData = _anaIntCache[code] ;
357  cacheData.S22bar.ResizeTo(S22bar) ;
358  cacheData.S22bar=S22bar ;
359  cacheData.S22det= S22.Determinant() ;
360  cacheData.pmap = map1 ;
361  cacheData.nint = map2.size() ;
362 
363  return cacheData ;
364 }
365 
366 
367 
368 ////////////////////////////////////////////////////////////////////////////////
369 /// Special case: generate all observables
370 
371 Int_t RooMultiVarGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, Bool_t /*staticInitOK*/) const
372 {
373  if (directVars.getSize()==_x.getSize()) {
374  generateVars.add(directVars) ;
375  return -1 ;
376  }
377 
378  Int_t nx = _x.getSize() ;
379  if (nx>127) {
380  // Warn that analytical integration is only provided for the first 127 observables
381  coutW(Integration) << "RooMultiVarGaussian::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << _x.getSize()
382  << " observables, partial internal generation is only implemented for the first 127 observables" << endl ;
383  nx=127 ;
384  }
385 
386  // Advertise partial generation over all permutations of observables
387  Int_t code(0) ;
388  BitBlock bits ;
389  for (int i=0 ; i<_x.getSize() ; i++) {
390  RooAbsArg* arg = directVars.find(_x.at(i)->GetName()) ;
391  if (arg) {
392  bits.setBit(i) ;
393 // code |= (1<<i) ;
394  generateVars.add(*arg) ;
395  }
396  }
397 
398  // Map BitBlock into return code
399  for (UInt_t i=0 ; i<_aicMap.size() ; i++) {
400  if (_aicMap[i]==bits) {
401  code = i+1 ;
402  }
403  }
404  if (code==0) {
405  _aicMap.push_back(bits) ;
406  code = _aicMap.size() ;
407  }
408 
409 
410  return code ;
411 }
412 
413 
414 
415 ////////////////////////////////////////////////////////////////////////////////
416 /// Clear the GenData cache as its content is not invariant under changes in
417 /// the mu vector.
418 
420 {
421  _genCache.clear() ;
422 
423 }
424 
425 
426 
427 
428 ////////////////////////////////////////////////////////////////////////////////
429 /// Retrieve generator config from cache
430 
432 {
433  GenData& gd = genData(code) ;
434  TMatrixD& TU = gd.UT ;
435  Int_t nobs = TU.GetNcols() ;
436  vector<int>& omap = gd.omap ;
437 
438  while(1) {
439 
440  // Create unit Gaussian vector
441  TVectorD xgen(nobs);
442  for(Int_t k= 0; k <nobs; k++) {
443  xgen(k)= RooRandom::gaussian();
444  }
445 
446  // Apply transformation matrix
447  xgen *= TU ;
448 
449  // Apply shift
450  if (code == -1) {
451 
452  // Simple shift if we generate all observables
453  xgen += gd.mu1 ;
454 
455  } else {
456 
457  // Non-generated observable dependent shift for partial generations
458 
459  // mubar = mu1 + S12 S22Inv ( x2 - mu2)
460  TVectorD mubar(gd.mu1) ;
461  TVectorD x2(gd.pmap.size()) ;
462  for (UInt_t i=0 ; i<gd.pmap.size() ; i++) {
463  x2(i) = ((RooAbsReal*)_x.at(gd.pmap[i]))->getVal() ;
464  }
465  mubar += gd.S12S22I * (x2 - gd.mu2) ;
466 
467  xgen += mubar ;
468 
469  }
470 
471  // Transfer values and check if values are in range
472  Bool_t ok(kTRUE) ;
473  for (int i=0 ; i<nobs ; i++) {
474  RooRealVar* xi = (RooRealVar*)_x.at(omap[i]) ;
475  if (xgen(i)<xi->getMin() || xgen(i)>xi->getMax()) {
476  ok = kFALSE ;
477  break ;
478  } else {
479  xi->setVal(xgen(i)) ;
480  }
481  }
482 
483  // If all values are in range, accept event and return
484  // otherwise retry
485  if (ok) {
486  break ;
487  }
488  }
489 
490  return;
491 }
492 
493 
494 
495 ////////////////////////////////////////////////////////////////////////////////
496 /// WVE -- CHECK THAT GENDATA IS VALID GIVEN CURRENT VALUES OF _MU
497 
499 {
500  // Check if cache entry was previously created
501  map<int,GenData>::iterator iter = _genCache.find(code) ;
502  if (iter != _genCache.end()) {
503  return iter->second ;
504  }
505 
506  // Create new entry
507  GenData& cacheData = _genCache[code] ;
508 
509  if (code==-1) {
510 
511  // Do eigen value decomposition
512  TDecompChol tdc(_cov) ;
513  tdc.Decompose() ;
514  TMatrixD U = tdc.GetU() ;
516 
517  // Fill cache data
518  cacheData.UT.ResizeTo(TU) ;
519  cacheData.UT = TU ;
520  cacheData.omap.resize(_x.getSize()) ;
521  for (int i=0 ; i<_x.getSize() ; i++) {
522  cacheData.omap[i] = i ;
523  }
524  syncMuVec() ;
525  cacheData.mu1.ResizeTo(_muVec) ;
526  cacheData.mu1 = _muVec ;
527 
528  } else {
529 
530  // Construct observables: map1 = generated, map2 = given
531  vector<int> map1, map2 ;
532  decodeCode(code,map2,map1) ;
533 
534  // Do block decomposition of covariance matrix
535  TMatrixDSym S11, S22 ;
536  TMatrixD S12, S21 ;
537  blockDecompose(_cov,map1,map2,S11,S12,S21,S22) ;
538 
539  // Constructed conditional matrix form
540  // -1
541  // F(X1|X2) --> CovI --> S22bar = S11 - S12 S22 S21
542  // -1
543  // --> mu --> mubar = mu1 + S12 S22 ( x2 - mu2)
544 
545  // Do eigenvalue decomposition
546  TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
547  TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
548 
549  // Do eigen value decomposition of S22bar
550  TDecompChol tdc(S22bar) ;
551  tdc.Decompose() ;
552  TMatrixD U = tdc.GetU() ;
554 
555  // Split mu vector into mu1 and mu2
556  TVectorD mu1(map1.size()),mu2(map2.size()) ;
557  syncMuVec() ;
558  for (UInt_t i=0 ; i<map1.size() ; i++) {
559  mu1(i) = _muVec(map1[i]) ;
560  }
561  for (UInt_t i=0 ; i<map2.size() ; i++) {
562  mu2(i) = _muVec(map2[i]) ;
563  }
564 
565  // Calculate rotation matrix for mu vector
566  TMatrixD S12S22Inv = S12 * S22Inv ;
567 
568  // Fill cache data
569  cacheData.UT.ResizeTo(TU) ;
570  cacheData.UT = TU ;
571  cacheData.omap = map1 ;
572  cacheData.pmap = map2 ;
573  cacheData.mu1.ResizeTo(mu1) ;
574  cacheData.mu2.ResizeTo(mu2) ;
575  cacheData.mu1 = mu1 ;
576  cacheData.mu2 = mu2 ;
577  cacheData.S12S22I.ResizeTo(S12S22Inv) ;
578  cacheData.S12S22I = S12S22Inv ;
579 
580  }
581 
582 
583  return cacheData ;
584 }
585 
586 
587 
588 
589 ////////////////////////////////////////////////////////////////////////////////
590 /// Decode analytical integration/generation code into index map of integrated/generated (map2)
591 /// and non-integrated/generated observables (map1)
592 
593 void RooMultiVarGaussian::decodeCode(Int_t code, vector<int>& map1, vector<int>& map2) const
594 {
595  if (code<0 || code> (Int_t)_aicMap.size()) {
596  cout << "RooMultiVarGaussian::decodeCode(" << GetName() << ") ERROR don't have bit pattern for code " << code << endl ;
597  throw string("RooMultiVarGaussian::decodeCode() ERROR don't have bit pattern for code") ;
598  }
599 
600  BitBlock b = _aicMap[code-1] ;
601  map1.clear() ;
602  map2.clear() ;
603  for (int i=0 ; i<_x.getSize() ; i++) {
604  if (b.getBit(i)) {
605  map2.push_back(i) ;
606  } else {
607  map1.push_back(i) ;
608  }
609  }
610 }
611 
612 
613 ////////////////////////////////////////////////////////////////////////////////
614 /// Block decomposition of covI according to given maps of observables
615 
616 void RooMultiVarGaussian::blockDecompose(const TMatrixD& input, const vector<int>& map1, const vector<int>& map2, TMatrixDSym& S11, TMatrixD& S12, TMatrixD& S21, TMatrixDSym& S22)
617 {
618  // Allocate and fill reordered covI matrix in 2x2 block structure
619 
620  S11.ResizeTo(map1.size(),map1.size()) ;
621  S12.ResizeTo(map1.size(),map2.size()) ;
622  S21.ResizeTo(map2.size(),map1.size()) ;
623  S22.ResizeTo(map2.size(),map2.size()) ;
624 
625  for (UInt_t i=0 ; i<map1.size() ; i++) {
626  for (UInt_t j=0 ; j<map1.size() ; j++)
627  S11(i,j) = input(map1[i],map1[j]) ;
628  for (UInt_t j=0 ; j<map2.size() ; j++)
629  S12(i,j) = input(map1[i],map2[j]) ;
630  }
631  for (UInt_t i=0 ; i<map2.size() ; i++) {
632  for (UInt_t j=0 ; j<map1.size() ; j++)
633  S21(i,j) = input(map2[i],map1[j]) ;
634  for (UInt_t j=0 ; j<map2.size() ; j++)
635  S22(i,j) = input(map2[i],map2[j]) ;
636  }
637 
638 }
639 
640 
642 {
643  if (ibit<32) { b0 |= (1<<ibit) ; return ; }
644  if (ibit<64) { b1 |= (1<<(ibit-32)) ; return ; }
645  if (ibit<96) { b2 |= (1<<(ibit-64)) ; return ; }
646  if (ibit<128) { b3 |= (1<<(ibit-96)) ; return ; }
647 }
648 
650 {
651  if (ibit<32) return (b0 & (1<<ibit)) ;
652  if (ibit<64) return (b1 & (1<<(ibit-32))) ;
653  if (ibit<96) return (b2 & (1<<(ibit-64))) ;
654  if (ibit<128) return (b3 & (1<<(ibit-96))) ;
655  return kFALSE ;
656 }
657 
659 {
660  if (b0 != other.b0) return kFALSE ;
661  if (b1 != other.b1) return kFALSE ;
662  if (b2 != other.b2) return kFALSE ;
663  if (b3 != other.b3) return kFALSE ;
664  return kTRUE ;
665 }
666 
667 
668 
669 
670 
const int nx
Definition: kalman.C:16
std::map< int, GenData > _genCache
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition: TVectorT.cxx:291
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Double_t evaluate() const
Do not persist.
const RooArgList & floatParsFinal() const
Definition: RooFitResult.h:109
static Double_t gaussian(TRandom *generator=randomGenerator())
Return a Gaussian random variable with mean 0 and variance 1.
Definition: RooRandom.cxx:110
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 double pi
#define cxcoutD(a)
Definition: RooMsgService.h:80
Int_t GetNrows() const
Definition: TVectorT.h:81
Bool_t contains(const RooAbsArg &var) const
virtual Double_t getMin(const char *name=0) const
Bool_t operator==(const BitBlock &other)
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:34
virtual TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1)
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
Definition: TMatrixT.cxx:1201
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, Bool_t staticInitOK=kTRUE) const
Special case: generate all observables.
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Handle full integral here.
std::vector< BitBlock > _aicMap
double sqrt(double)
std::map< int, AnaIntData > _anaIntCache
static const double x2[5]
Double_t x[n]
Definition: legend1.C:17
virtual Double_t Determinant() const
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
double pow(double, double)
void initGenerator(Int_t code)
Clear the GenData cache as its content is not invariant under changes in the mu vector.
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
TMatrixT< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant.
Definition: TMatrixT.cxx:1387
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:202
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
ClassImp(RooMultiVarGaussian)
RooAbsArg * find(const char *name) const
Find object with given name in list.
return
Definition: TBase64.cxx:62
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
void setConstant(Bool_t value=kTRUE)
const TMatrixD & GetU() const
Definition: TDecompChol.h:49
void generateEvent(Int_t code)
Retrieve generator config from cache.
unsigned int UInt_t
Definition: RtypesCore.h:42
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
AnaIntData & anaIntData(Int_t code) const
Check if cache entry was previously created.
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
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
virtual TObject * Clone(const char *newname=0) const
Make a clone of an object using the Streamer facility.
Definition: RooAbsArg.h:75
Int_t GetNcols() const
Definition: TMatrixTBase.h:137
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...
virtual Bool_t addOwned(RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::addOwned()
virtual Bool_t Decompose()
Matrix A is decomposed in component U so that A = U^T * U If the decomposition succeeds, bit kDecomposed is set , otherwise kSingular.
#define name(a, b)
Definition: linkTestLib0.cpp:5
TMatrixTSym< Element > & Invert(Double_t *det=0)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
virtual Bool_t remove(const RooAbsArg &var, Bool_t silent=kFALSE, Bool_t matchByNameOnly=kFALSE)
Remove the specified argument from our list.
virtual Double_t getMax(const char *name=0) const
RooConstVar & RooConst(Double_t val)
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
double exp(double)
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448
GenData & genData(Int_t code) const
WVE – CHECK THAT GENDATA IS VALID GIVEN CURRENT VALUES OF _MU.