Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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\file RooMultiVarGaussian.cxx
19\class RooMultiVarGaussian
20\ingroup Roofitcore
21
22Multivariate Gaussian p.d.f. with correlations
23**/
24
25#include "Riostream.h"
26#include <cmath>
27
28#include "RooMultiVarGaussian.h"
29#include "RooAbsReal.h"
30#include "RooRealVar.h"
31#include "RooRandom.h"
32#include "RooMath.h"
33#include "RooGlobalFunc.h"
34#include "RooConstVar.h"
35#include "TDecompChol.h"
36#include "RooFitResult.h"
37
38using std::string, std::list, std::map, std::vector, std::cout, std::endl;
39
41
42////////////////////////////////////////////////////////////////////////////////
43
44RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title,
45 const RooArgList& xvec, const RooArgList& mu, const TMatrixDBase& cov) :
46 RooAbsPdf(name,title),
47 _x("x","Observables",this,true,false),
48 _mu("mu","Offset vector",this,true,false),
49 _cov{cov.GetNrows()},
50 _covI{cov.GetNrows()},
51 _z(4)
52{
53 if(!cov.IsSymmetric()) {
54 std::stringstream errorMsg;
55 errorMsg << "RooMultiVarGaussian::RooMultiVarGaussian(" << GetName()
56 << ") input covariance matrix is not symmetric!";
57 coutE(InputArguments) << errorMsg.str() << std::endl;
58 throw std::invalid_argument(errorMsg.str().c_str());
59 }
60
61 _cov.SetSub(0, cov);
62 _covI.SetSub(0, cov);
63
64 _x.add(xvec) ;
65
66 _mu.add(mu) ;
67
69
70 // Invert covariance matrix
71 _covI.Invert() ;
72}
73
74
75////////////////////////////////////////////////////////////////////////////////
76
77RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title, const RooArgList &xvec,
78 const RooFitResult &fr, bool reduceToConditional)
79 : RooAbsPdf(name, title),
80 _x("x", "Observables", this, true, false),
81 _mu("mu", "Offset vector", this, true, false),
82 _cov(reduceToConditional ? fr.conditionalCovarianceMatrix(xvec) : fr.reducedCovarianceMatrix(xvec)),
83 _covI(_cov),
84 _det(_cov.Determinant()),
85 _z(4)
86{
87
88 // Fill mu vector with constant RooRealVars
89 list<string> munames ;
90 const RooArgList& fpf = fr.floatParsFinal() ;
91 for (std::size_t i=0 ; i<fpf.size() ; i++) {
92 if (xvec.find(fpf.at(i)->GetName())) {
93 std::unique_ptr<RooRealVar> parclone{static_cast<RooRealVar*>(fpf.at(i)->Clone(Form("%s_centralvalue",fpf.at(i)->GetName())))};
94 parclone->setConstant(true) ;
95 _mu.addOwned(std::move(parclone));
96 munames.push_back(fpf.at(i)->GetName()) ;
97 }
98 }
99
100 // Fill X vector in same order as mu vector
101 for (list<string>::iterator iter=munames.begin() ; iter!=munames.end() ; ++iter) {
102 RooRealVar* xvar = static_cast<RooRealVar*>(xvec.find(iter->c_str())) ;
103 _x.add(*xvar) ;
104 }
105
106 // Invert covariance matrix
107 _covI.Invert() ;
108
109}
110
111namespace {
112
113RooArgList muFromVector(TVectorD const &mu)
114{
115 RooArgList out;
116 for (int i = 0; i < mu.GetNrows(); i++) {
117 out.add(RooFit::RooConst(mu(i)));
118 }
119 return out;
120}
121
122RooArgList muConstants(std::size_t n)
123{
125 for (std::size_t i = 0; i < n; i++) {
126 out.add(RooFit::RooConst(0));
127 }
128 return out;
129}
130
131} // namespace
132
133
134////////////////////////////////////////////////////////////////////////////////
135
136RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title, const RooArgList &xvec,
137 const TVectorD &mu, const TMatrixDBase &cov)
138 : RooMultiVarGaussian{name, title, xvec, muFromVector(mu), cov}
139{
140}
141
142////////////////////////////////////////////////////////////////////////////////
143
144RooMultiVarGaussian::RooMultiVarGaussian(const char *name, const char *title, const RooArgList &xvec,
145 const TMatrixDBase &cov)
146 : RooMultiVarGaussian{name, title, xvec, muConstants(xvec.size()), cov}
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 (std::size_t i=0 ; i<_mu.size() ; i++) {
167 _muVec[i] = static_cast<RooAbsReal*>(_mu.at(i))->getVal() ;
168 }
169}
170
171
172////////////////////////////////////////////////////////////////////////////////
173/// Represent observables as vector
174
176{
177 TVectorD x(_x.size()) ;
178 for (std::size_t i=0 ; i<_x.size() ; i++) {
179 x[i] = static_cast<RooAbsReal*>(_x.at(i))->getVal() ;
180 }
181
182 // Calculate return value
183 syncMuVec() ;
184 TVectorD x_min_mu = x - _muVec ;
185
186 double alpha = x_min_mu * (_covI * x_min_mu) ;
187 return exp(-0.5*alpha) ;
188}
189
190
191
192////////////////////////////////////////////////////////////////////////////////
193
194Int_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 (std::size_t i=0 ; i<_x.size() ; i++) {
200 if (allVars.contains(*_x.at(i))) {
201 allVars.remove(*_mu.at(i),true,true) ;
202 }
203 }
204
205
206 // Analytical integral known over all observables
207 if (allVars.size()==_x.size() && !rangeName) {
208 analVars.add(allVars) ;
209 return -1 ;
210 }
211
212 Int_t code(0) ;
213
214 Int_t nx = _x.size() ;
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.size()
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 anyBits(false) ;
226 syncMuVec() ;
227 for (std::size_t i=0 ; i<_x.size() ; 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 = static_cast<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 = true ;
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 = static_cast<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 = true ;
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
288double RooMultiVarGaussian::analyticalIntegral(Int_t code, const char* /*rangeName*/) const
289{
290 if (code==-1) {
291 return pow(2*3.14159268,_x.size()/2.)*sqrt(std::abs(_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) = (static_cast<RooAbsReal*>(_x.at(aid.pmap[i])))->getVal() - _muVec(aid.pmap[i]) ;
304 }
305
306 // Calculate partial integral
307 double ret = pow(2*3.14159268,aid.nint/2.)/sqrt(std::abs(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;
328 vector<int> map2;
329 decodeCode(code,map1,map2) ;
330
331 // Rearrage observables so that all non-integrated observables
332 // go first (preserving relative order) and all integrated observables
333 // go last (preserving relative order)
334 TMatrixDSym S11;
335 TMatrixDSym S22;
336 TMatrixD S12;
337 TMatrixD S21;
338 blockDecompose(_covI,map1,map2,S11,S12,S21,S22) ;
339
340 // Begin calculation of partial integrals
341 // ___
342 // sqrt(2pi)^(#intObs) (-0.5 * u1T S22 u1 )
343 // I = ------------------- * e
344 // sqrt(|det(S22)|)
345 // ___
346 // Where S22 is the sub-matrix of covI for the integrated observables and S22
347 // is the Schur complement of S22
348 // ___ -1
349 // S22 = S11 - S12 * S22 * S21
350 //
351 // and u1 is the vector of non-integrated observables
352
353 // Calculate Schur complement S22bar
354 TMatrixD S22inv(S22) ;
355 S22inv.Invert() ;
356 TMatrixD S22bar = S11 - S12*S22inv*S21 ;
357
358 // Create new cache entry
359 AnaIntData& cacheData = _anaIntCache[code] ;
360 cacheData.S22bar.ResizeTo(S22bar) ;
361 cacheData.S22bar=S22bar ;
362 cacheData.S22det= S22.Determinant() ;
363 cacheData.pmap = map1 ;
364 cacheData.nint = map2.size() ;
365
366 return cacheData ;
367}
368
369
370
371////////////////////////////////////////////////////////////////////////////////
372/// Special case: generate all observables
373
374Int_t RooMultiVarGaussian::getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool /*staticInitOK*/) const
375{
376 if (directVars.size()==_x.size()) {
377 generateVars.add(directVars) ;
378 return -1 ;
379 }
380
381 Int_t nx = _x.size() ;
382 if (nx>127) {
383 // Warn that analytical integration is only provided for the first 127 observables
384 coutW(Integration) << "RooMultiVarGaussian::getGenerator(" << GetName() << ") WARNING: p.d.f. has " << _x.size()
385 << " observables, partial internal generation is only implemented for the first 127 observables" << endl ;
386 nx=127 ;
387 }
388
389 // Advertise partial generation over all permutations of observables
390 Int_t code(0) ;
391 BitBlock bits ;
392 for (std::size_t i=0 ; i<_x.size() ; i++) {
393 RooAbsArg* arg = directVars.find(_x.at(i)->GetName()) ;
394 if (arg) {
395 bits.setBit(i) ;
396// code |= (1<<i) ;
397 generateVars.add(*arg) ;
398 }
399 }
400
401 // Map BitBlock into return code
402 for (UInt_t i=0 ; i<_aicMap.size() ; i++) {
403 if (_aicMap[i]==bits) {
404 code = i+1 ;
405 }
406 }
407 if (code==0) {
408 _aicMap.push_back(bits) ;
409 code = _aicMap.size() ;
410 }
411
412
413 return code ;
414}
415
416
417
418////////////////////////////////////////////////////////////////////////////////
419/// Clear the GenData cache as its content is not invariant under changes in
420/// the mu vector.
421
423{
424 _genCache.clear() ;
425
426}
427
428
429
430
431////////////////////////////////////////////////////////////////////////////////
432/// Retrieve generator config from cache
433
435{
436 GenData& gd = genData(code) ;
437 TMatrixD& TU = gd.UT ;
438 Int_t nobs = TU.GetNcols() ;
439 vector<int>& omap = gd.omap ;
440
441 while(true) {
442
443 // Create unit Gaussian vector
444 TVectorD xgen(nobs);
445 for(Int_t k= 0; k <nobs; k++) {
446 xgen(k)= RooRandom::gaussian();
447 }
448
449 // Apply transformation matrix
450 xgen *= TU ;
451
452 // Apply shift
453 if (code == -1) {
454
455 // Simple shift if we generate all observables
456 xgen += gd.mu1 ;
457
458 } else {
459
460 // Non-generated observable dependent shift for partial generations
461
462 // mubar = mu1 + S12 S22Inv ( x2 - mu2)
463 TVectorD mubar(gd.mu1) ;
464 TVectorD x2(gd.pmap.size()) ;
465 for (UInt_t i=0 ; i<gd.pmap.size() ; i++) {
466 x2(i) = (static_cast<RooAbsReal*>(_x.at(gd.pmap[i])))->getVal() ;
467 }
468 mubar += gd.S12S22I * (x2 - gd.mu2) ;
469
470 xgen += mubar ;
471
472 }
473
474 // Transfer values and check if values are in range
475 bool ok(true) ;
476 for (int i=0 ; i<nobs ; i++) {
477 RooRealVar* xi = static_cast<RooRealVar*>(_x.at(omap[i])) ;
478 if (xgen(i)<xi->getMin() || xgen(i)>xi->getMax()) {
479 ok = false ;
480 break ;
481 } else {
482 xi->setVal(xgen(i)) ;
483 }
484 }
485
486 // If all values are in range, accept event and return
487 // otherwise retry
488 if (ok) {
489 break ;
490 }
491 }
492
493 return;
494}
495
496
497
498////////////////////////////////////////////////////////////////////////////////
499/// WVE -- CHECK THAT GENDATA IS VALID GIVEN CURRENT VALUES OF _MU
500
502{
503 // Check if cache entry was previously created
504 map<int,GenData>::iterator iter = _genCache.find(code) ;
505 if (iter != _genCache.end()) {
506 return iter->second ;
507 }
508
509 // Create new entry
510 GenData& cacheData = _genCache[code] ;
511
512 if (code==-1) {
513
514 // Do eigen value decomposition
515 TDecompChol tdc(_cov) ;
516 tdc.Decompose() ;
517 TMatrixD U = tdc.GetU() ;
519
520 // Fill cache data
521 cacheData.UT.ResizeTo(TU) ;
522 cacheData.UT = TU ;
523 cacheData.omap.resize(_x.size()) ;
524 for (std::size_t i=0 ; i<_x.size() ; i++) {
525 cacheData.omap[i] = i ;
526 }
527 syncMuVec() ;
528 cacheData.mu1.ResizeTo(_muVec) ;
529 cacheData.mu1 = _muVec ;
530
531 } else {
532
533 // Construct observables: map1 = generated, map2 = given
534 vector<int> map1;
535 vector<int> map2;
536 decodeCode(code,map2,map1) ;
537
538 // Do block decomposition of covariance matrix
539 TMatrixDSym S11;
540 TMatrixDSym S22;
541 TMatrixD S12;
542 TMatrixD S21;
543 blockDecompose(_cov,map1,map2,S11,S12,S21,S22) ;
544
545 // Constructed conditional matrix form
546 // -1
547 // F(X1|X2) --> CovI --> S22bar = S11 - S12 S22 S21
548 // -1
549 // --> mu --> mubar = mu1 + S12 S22 ( x2 - mu2)
550
551 // Do eigenvalue decomposition
552 TMatrixD S22Inv(TMatrixD::kInverted,S22) ;
553 TMatrixD S22bar = S11 - S12 * (S22Inv * S21) ;
554
555 // Do eigen value decomposition of S22bar
556 TDecompChol tdc(S22bar) ;
557 tdc.Decompose() ;
558 TMatrixD U = tdc.GetU() ;
560
561 // Split mu vector into mu1 and mu2
562 TVectorD mu1(map1.size());
563 TVectorD mu2(map2.size());
564 syncMuVec() ;
565 for (UInt_t i=0 ; i<map1.size() ; i++) {
566 mu1(i) = _muVec(map1[i]) ;
567 }
568 for (UInt_t i=0 ; i<map2.size() ; i++) {
569 mu2(i) = _muVec(map2[i]) ;
570 }
571
572 // Calculate rotation matrix for mu vector
573 TMatrixD S12S22Inv = S12 * S22Inv ;
574
575 // Fill cache data
576 cacheData.UT.ResizeTo(TU) ;
577 cacheData.UT = TU ;
578 cacheData.omap = map1 ;
579 cacheData.pmap = map2 ;
580 cacheData.mu1.ResizeTo(mu1) ;
581 cacheData.mu2.ResizeTo(mu2) ;
582 cacheData.mu1 = mu1 ;
583 cacheData.mu2 = mu2 ;
584 cacheData.S12S22I.ResizeTo(S12S22Inv) ;
585 cacheData.S12S22I = S12S22Inv ;
586
587 }
588
589
590 return cacheData ;
591}
592
593
594
595
596////////////////////////////////////////////////////////////////////////////////
597/// Decode analytical integration/generation code into index map of integrated/generated (map2)
598/// and non-integrated/generated observables (map1)
599
600void RooMultiVarGaussian::decodeCode(Int_t code, vector<int>& map1, vector<int>& map2) const
601{
602 if (code<0 || code> (Int_t)_aicMap.size()) {
603 cout << "RooMultiVarGaussian::decodeCode(" << GetName() << ") ERROR don't have bit pattern for code " << code << endl ;
604 throw string("RooMultiVarGaussian::decodeCode() ERROR don't have bit pattern for code") ;
605 }
606
607 BitBlock b = _aicMap[code-1] ;
608 map1.clear() ;
609 map2.clear() ;
610 for (std::size_t i=0 ; i<_x.size() ; i++) {
611 if (b.getBit(i)) {
612 map2.push_back(i) ;
613 } else {
614 map1.push_back(i) ;
615 }
616 }
617}
618
619
620////////////////////////////////////////////////////////////////////////////////
621/// Block decomposition of covI according to given maps of observables
622
623void RooMultiVarGaussian::blockDecompose(const TMatrixD& input, const vector<int>& map1, const vector<int>& map2, TMatrixDSym& S11, TMatrixD& S12, TMatrixD& S21, TMatrixDSym& S22)
624{
625 // Allocate and fill reordered covI matrix in 2x2 block structure
626
627 S11.ResizeTo(map1.size(),map1.size()) ;
628 S12.ResizeTo(map1.size(),map2.size()) ;
629 S21.ResizeTo(map2.size(),map1.size()) ;
630 S22.ResizeTo(map2.size(),map2.size()) ;
631
632 for (UInt_t i=0 ; i<map1.size() ; i++) {
633 for (UInt_t j=0 ; j<map1.size() ; j++)
634 S11(i,j) = input(map1[i],map1[j]) ;
635 for (UInt_t j=0 ; j<map2.size() ; j++)
636 S12(i,j) = input(map1[i],map2[j]) ;
637 }
638 for (UInt_t i=0 ; i<map2.size() ; i++) {
639 for (UInt_t j=0 ; j<map1.size() ; j++)
640 S21(i,j) = input(map2[i],map1[j]) ;
641 for (UInt_t j=0 ; j<map2.size() ; j++)
642 S22(i,j) = input(map2[i],map2[j]) ;
643 }
644
645}
646
647
649{
650 if (ibit<32) { b0 |= (1<<ibit) ; return ; }
651 if (ibit<64) { b1 |= (1<<(ibit-32)) ; return ; }
652 if (ibit<96) { b2 |= (1<<(ibit-64)) ; return ; }
653 if (ibit<128) { b3 |= (1<<(ibit-96)) ; return ; }
654}
655
657{
658 if (ibit<32) return (b0 & (1<<ibit)) ;
659 if (ibit<64) return (b1 & (1<<(ibit-32))) ;
660 if (ibit<96) return (b2 & (1<<(ibit-64))) ;
661 if (ibit<128) return (b3 & (1<<(ibit-96))) ;
662 return false ;
663}
664
666{
667 if (lhs.b0 != rhs.b0)
668 return false;
669 if (lhs.b1 != rhs.b1)
670 return false;
671 if (lhs.b2 != rhs.b2)
672 return false;
673 if (lhs.b3 != rhs.b3)
674 return false;
675 return true;
676}
#define b(i)
Definition RSha256.hxx:100
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define cxcoutD(a)
#define coutW(a)
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:377
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void input
Option_t Option_t TPoint TPoint const char x2
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
TObject * Clone(const char *newname=nullptr) const override
Make a clone of an object using the Streamer facility.
Definition RooAbsArg.h:89
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
bool contains(const RooAbsArg &var) const
Check if collection contains an argument with the same name as var.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
void setConstant(bool value=true)
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
bool operator==(double value) const
Equality operator comparing to a double.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
bool addOwned(RooAbsArg &var, bool silent=false) override
Overloaded RooCollection_t::addOwned() method insert object into owning set and registers object as s...
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
RooFitResult is a container class to hold the input and output of a PDF fit to a dataset.
const RooArgList & floatParsFinal() const
Return list of floating parameters after fit.
Multivariate Gaussian p.d.f.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Handle full integral here.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
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...
std::vector< BitBlock > _aicMap
!
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.
AnaIntData & anaIntData(Int_t code) const
Check if cache entry was previously created.
double evaluate() const override
Do not persist.
std::map< int, GenData > _genCache
!
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Special case: generate all observables.
void initGenerator(Int_t code) override
Clear the GenData cache as its content is not invariant under changes in the mu vector.
GenData & genData(Int_t code) const
WVE – CHECK THAT GENDATA IS VALID GIVEN CURRENT VALUES OF _MU.
std::map< int, AnaIntData > _anaIntCache
!
void generateEvent(Int_t code) override
Retrieve generator config from cache.
static double gaussian(TRandom *generator=randomGenerator())
Return a Gaussian random variable with mean 0 and variance 1.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
Cholesky Decomposition class.
Definition TDecompChol.h:25
Bool_t Decompose() override
Matrix A is decomposed in component U so that A = U^T * U If the decomposition succeeds,...
const TMatrixD & GetU() const
Definition TDecompChol.h:45
Int_t GetNcols() const
virtual Bool_t IsSymmetric() const
Check whether matrix is symmetric.
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
TMatrixTSym< Element > & SetSub(Int_t row_lwb, const TMatrixTBase< Element > &source)
Insert matrix source starting at [row_lwb][row_lwb], thereby overwriting the part [row_lwb....
Double_t Determinant() const override
TMatrixTSym< Element > & Invert(Double_t *det=nullptr)
Invert the matrix and calculate its determinant Notice that the LU decomposition is used instead of B...
TMatrixTBase< Element > & ResizeTo(Int_t nrows, Int_t ncols, Int_t=-1) override
Set size of the matrix to nrows x ncols New dynamic elements are created, the overlapping part of the...
TMatrixT< Element > & Invert(Double_t *det=nullptr)
Invert the matrix and calculate its determinant.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
TVectorT< Element > & ResizeTo(Int_t lwb, Int_t upb)
Resize the vector to [lwb:upb] .
Definition TVectorT.cxx:294
Int_t GetNrows() const
Definition TVectorT.h:73
RooConstVar & RooConst(double val)
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16