Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBarlowBeestonLL.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id$
2// Author: Kyle Cranmer, George Lewis
3/*************************************************************************
4 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
5 * All rights reserved. *
6 * *
7 * For the licensing terms see $ROOTSYS/LICENSE. *
8 * For the list of contributors see $ROOTSYS/README/CREDITS. *
9 *************************************************************************/
10
11//////////////////////////////////////////////////////////////////////////////
12/** \class RooStats::HistFactory::RooBarlowBeestonLL
13 * \ingroup HistFactory
14//
15// Class RooBarlowBeestonLL implements the profile likelihood estimator for
16// a given likelihood and set of parameters of interest. The value return by
17// RooBarlowBeestonLL is the input likelihood nll minimized w.r.t all nuisance parameters
18// (which are all parameters except for those listed in the constructor) minus
19// the -log(L) of the best fit. Note that this function is slow to evaluate
20// as a MIGRAD minimization step is executed for each function evaluation
21*/
22
23#include <stdexcept>
24#include <cmath>
25#include <iostream>
26
28#include "RooAbsReal.h"
29#include "RooMsgService.h"
30#include "RooRealVar.h"
31
33#include "RooCategory.h"
34#include "RooSimultaneous.h"
35#include "RooArgList.h"
36
39
40#include <TMath.h>
41
42using std::runtime_error;
43
45
46
47////////////////////////////////////////////////////////////////////////////////
48
49/// Default constructor. Should only be used by proof.
51{
52}
53
54
55////////////////////////////////////////////////////////////////////////////////
56
58 RooAbsReal& nllIn /*, const RooArgSet& observables*/) :
59 RooAbsReal(name,title),
60 _nll("input","-log(L) function",this,nllIn)
61 // _obs("paramOfInterest","Parameters of interest",this),
62 // _par("nuisanceParam","Nuisance parameters",this,false,false),
63{
64 // Constructor of profile likelihood given input likelihood nll w.r.t
65 // the given set of variables. The input log likelihood is minimized w.r.t
66 // to all other variables of the likelihood at each evaluation and the
67 // value of the global log likelihood minimum is always subtracted.
68
69 // Determine actual parameters and observables
70 /*
71 std::unique_ptr<RooArgSet> actualObs{nllIn.getObservables(observables)};
72 std::unique_ptr<RooArgSet> actualPars{nllIn.getParameters(observables)};
73
74 _obs.add(*actualObs) ;
75 _par.add(*actualPars) ;
76 */
77}
78
79
80
81////////////////////////////////////////////////////////////////////////////////
82
84 RooAbsReal(other,name),
85 _nll("nll",this,other._nll),
86 // _obs("obs",this,other._obs),
87 // _par("par",this,other._par),
88 _paramFixed(other._paramFixed)
89{
90 // Copy constructor
91
92 // _paramAbsMin.addClone(other._paramAbsMin) ;
93 // _obsAbsMin.addClone(other._obsAbsMin) ;
94
95}
96
97
98////////////////////////////////////////////////////////////////////////////////
99
101 for (auto const *var : static_range_cast<RooRealVar *>(*bin_center)) {
102 RooRealVar* target = static_cast<RooRealVar*>(observables->find(var->GetName())) ;
103 target->setVal(var->getVal()) ;
104 }
105}
106
107
108////////////////////////////////////////////////////////////////////////////////
109
111 bool verbose=false;
112
113 if(!_data) {
114 std::cout << "Error: Must initialize data before initializing cache" << std::endl;
115 throw std::runtime_error("Uninitialized Data");
116 }
117 if(!_pdf) {
118 std::cout << "Error: Must initialize model pdf before initializing cache" << std::endl;
119 throw std::runtime_error("Uninitialized model pdf");
120 }
121
122 // Get the data bin values for all channels and bins
123 std::map< std::string, std::vector<double> > ChannelBinDataMap;
124 getDataValuesForObservables( ChannelBinDataMap, _data, _pdf );
125
126 // Get a list of constraint terms
127 RooArgList obsTerms;
128 RooArgList constraints;
129 RooArgSet* obsSet = std::unique_ptr<RooArgSet>{_pdf->getObservables(*_data)}.release();
130 FactorizeHistFactoryPdf(*obsSet, *_pdf, obsTerms, constraints);
131
132 if( obsTerms.empty() ) {
133 std::cout << "Error: Found no observable terms with pdf: " << _pdf->GetName()
134 << " using dataset: " << _data->GetName() << std::endl;
135 return;
136 }
137 if( constraints.empty() ) {
138 std::cout << "Error: Found no constraint terms with pdf: " << _pdf->GetName()
139 << " using dataset: " << _data->GetName() << std::endl;
140 return;
141 }
142
143 // Loop over the channels
144 auto simPdf = static_cast<RooSimultaneous*>(_pdf);
145 auto channelCat = static_cast<RooCategory const*>(&simPdf->indexCat());
146 for (const auto& nameIdx : *channelCat) {
147
148 // Warning: channel cat name is not necessarily the same name
149 // as the pdf's (for example, if someone does edits)
150 RooAbsPdf* channelPdf = simPdf->getPdf(nameIdx.first.c_str());
151 std::string channel_name = channelPdf->GetName();
152
153 // First, we check if this channel uses Stat Uncertainties:
154 RooArgList* gammas = nullptr;
155 ParamHistFunc* param_func=nullptr;
156 bool hasStatUncert = getStatUncertaintyFromChannel( channelPdf, param_func, gammas );
157 if( ! hasStatUncert ) {
158 if(verbose) {
159 std::cout << "Channel: " << channel_name
160 << " doesn't have statistical uncertainties"
161 << std::endl;
162 }
163 continue;
164 }
165 else {
166 if(verbose) std::cout << "Found ParamHistFunc: " << param_func->GetName() << std::endl;
167 }
168
169 // Now, loop over the bins in this channel
170 // To Do: Check that the index convention
171 // still works for 2-d (ie matches the
172 // convention in ParamHistFunc, etc)
173 int num_bins = param_func->numBins();
174
175 // Initialize the vector to the number of bins, allowing
176 // us to skip gamma's if they are constant
177
178 std::vector<BarlowCache> temp_cache( num_bins );
179 bool channel_has_stat_uncertainty=false;
180
181 for( Int_t bin_index = 0; bin_index < num_bins; ++bin_index ) {
182
183 // Create a cache object
184 BarlowCache cache;
185
186 // Get the gamma for this bin, and skip if it's constant
187 RooRealVar* gamma_stat = dynamic_cast<RooRealVar*>(&(param_func->getParameter(bin_index)));
188 if(!gamma_stat) throw std::runtime_error("ParamHistFunc contains non-RooRealVar, not supported in RooBarlowBeestonLL");
189 if( gamma_stat->isConstant() ) {
190 if(verbose) std::cout << "Ignoring constant gamma: " << gamma_stat->GetName() << std::endl;
191 continue;
192 }
193 else {
194 cache.hasStatUncert=true;
195 channel_has_stat_uncertainty=true;
196 cache.gamma = gamma_stat;
197 _statUncertParams.insert( gamma_stat->GetName() );
198 }
199
200 // Store a snapshot of the bin center
201 RooArgSet* bin_center = (RooArgSet*) param_func->get( bin_index )->snapshot();
202 cache.bin_center = bin_center;
203 cache.observables = obsSet;
204
205 cache.binVolume = param_func->binVolume();
206
207 // Delete this part, simply a comment
208 RooArgList obs_list( *(cache.bin_center) );
209
210 // Get the gamma's constraint term
211 RooAbsReal* pois_mean = nullptr;
212 RooRealVar* tau = nullptr;
213 getStatUncertaintyConstraintTerm( &constraints, gamma_stat, pois_mean, tau );
214 if( !tau || !pois_mean ) {
215 std::cout << "Failed to find pois mean or tau parameter for " << gamma_stat->GetName() << std::endl;
216 }
217 else {
218 if (verbose) {
219 std::cout << "Found pois mean and tau for parameter: " << gamma_stat->GetName() << " tau: " << tau->GetName()
220 << " " << tau->getVal() << " pois_mean: " << pois_mean->GetName() << " " << pois_mean->getVal()
221 << std::endl;
222 }
223 }
224
225 cache.tau = tau;
226 cache.nom_pois_mean = pois_mean;
227
228 // Get the RooRealSumPdf
229 RooAbsPdf* sum_pdf = getSumPdfFromChannel( channelPdf );
230 if( sum_pdf == nullptr ) {
231 std::cout << "Failed to find RooRealSumPdf in channel " << channel_name
232 << ", therefor skipping this channel for analytic uncertainty minimization"
233 << std::endl;
234 channel_has_stat_uncertainty=false;
235 break;
236 }
237 cache.sumPdf = sum_pdf;
238
239 // And set the data value for this bin
240 if( ChannelBinDataMap.find(channel_name) == ChannelBinDataMap.end() ) {
241 std::cout << "Error: channel with name: " << channel_name
242 << " not found in BinDataMap" << std::endl;
243 throw runtime_error("BinDataMap");
244 }
245 double nData = ChannelBinDataMap[channel_name].at(bin_index);
246 cache.nData = nData;
247
248 temp_cache.at(bin_index) = cache;
249 //_barlowCache[channel_name].at(bin_index) = cache;
250
251 } // End: Loop over bins
252
253 if( channel_has_stat_uncertainty ) {
254 std::cout << "Adding channel: " << channel_name
255 << " to the barlow cache" << std::endl;
256 _barlowCache[channel_name] = temp_cache;
257 }
258
259
260 } // End: Loop over channels
261
262
263
264 // Successfully initialized the cache
265 // Printing some info
266 /*
267 std::map< std::string, std::vector< BarlowCache > >::iterator iter_cache;
268 for( iter_cache = _barlowCache.begin(); iter_cache != _barlowCache.end(); ++iter_cache ) {
269
270 std::string channel_name = (*iter_cache).first;
271 std::vector< BarlowCache >& channel_cache = (*iter_cache).second;
272
273
274 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
275 BarlowCache& bin_cache = channel_cache.at(i);
276
277
278 RooRealVar* gamma = bin_cache.gamma;
279 RooRealVar* tau = bin_cache.tau;
280 RooAbsReal* pois_mean = bin_cache.nom_pois_mean;
281 RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
282 double binVolume = bin_cache.binVolume;
283
284
285 if( !bin_cache.hasStatUncert ) {
286 std::cout << "Barlow Cache for Channel: " << channel_name
287 << " Bin: " << i
288 << " NO STAT UNCERT"
289 << std::endl;
290 }
291 else {
292 std::cout << "Barlow Cache for Channel: " << channel_name
293 << " Bin: " << i
294 << " gamma: " << gamma->GetName()
295 << " tau: " << tau->GetName()
296 << " pois_mean: " << pois_mean->GetName()
297 << " sum_pdf: " << sum_pdf->GetName()
298 << " binVolume: " << binVolume
299 << std::endl;
300 }
301
302 }
303 }
304 */
305
306}
307
308
309////////////////////////////////////////////////////////////////////////////////
310
312 RooArgSet& outputSet,
313 bool stripDisconnected) const {
314 bool errorInBaseCall = RooAbsArg::getParameters( depList, outputSet, stripDisconnected );
315
316 RooArgSet toRemove;
317 toRemove.reserve( _statUncertParams.size());
318
319 for (auto const& arg : outputSet) {
320
321 // If there is a gamma in the name,
322 // strip it from the list of dependencies
323
324 if( _statUncertParams.find(arg->GetName()) != _statUncertParams.end() ) {
325 toRemove.add( *arg );
326 }
327 }
328
329 for( auto& arg : toRemove) outputSet.remove( *arg, true );
330
331 return errorInBaseCall || false;
332
333}
334
335
336/*
337////////////////////////////////////////////////////////////////////////////////
338
339const RooArgSet& RooStats::HistFactory::RooBarlowBeestonLL::bestFitParams() const
340{
341 validateAbsMin() ;
342 return _paramAbsMin ;
343}
344
345
346////////////////////////////////////////////////////////////////////////////////
347
348const RooArgSet& RooStats::HistFactory::RooBarlowBeestonLL::bestFitObs() const
349{
350 validateAbsMin() ;
351 return _obsAbsMin ;
352}
353*/
354
355
356
357////////////////////////////////////////////////////////////////////////////////
358/// Optimized implementation of createProfile for profile likelihoods.
359/// Return profile of original function in terms of stated parameters
360/// of interest rather than profiling recursively.
361
362/*
363RooAbsReal* RooStats::HistFactory::RooBarlowBeestonLL::createProfile(const RooArgSet& paramsOfInterest)
364{
365 return nll().createProfile(paramsOfInterest) ;
366}
367*/
368
369
370/*
371void RooStats::HistFactory::RooBarlowBeestonLL::FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) const {
372 // utility function to factorize constraint terms from a pdf
373 // (from G. Petrucciani)
374 const std::type_info & id = typeid(pdf);
375 if (id == typeid(RooProdPdf)) {
376 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
377 RooArgList list(prod->pdfList());
378 for (int i = 0, n = list.size(); i < n; ++i) {
379 RooAbsPdf *pdfi = (RooAbsPdf *) list.at(i);
380 FactorizePdf(observables, *pdfi, obsTerms, constraints);
381 }
382 } else if (id == typeid(RooSimultaneous) ) { //|| id == typeid(RooSimultaneousOpt)) {
383 RooSimultaneous *sim = dynamic_cast<RooSimultaneous *>(&pdf);
384 std::unique_ptr<RooAbsCategoryLValue> cat{(RooAbsCategoryLValue *) sim->indexCat().Clone()};
385 for (int ic = 0, nc = cat->numBins((const char *)0); ic < nc; ++ic) {
386 cat->setBin(ic);
387 FactorizePdf(observables, *sim->getPdf(cat->getLabel()), obsTerms, constraints);
388 }
389 } else if (pdf.dependsOn(observables)) {
390 if (!obsTerms.contains(pdf)) obsTerms.add(pdf);
391 } else {
392 if (!constraints.contains(pdf)) constraints.add(pdf);
393 }
394}
395*/
396
397
398
399////////////////////////////////////////////////////////////////////////////////
400
402{
403 // Loop over the channels (keys to the map)
404 //clock_t time_before_setVal, time_after_setVal;
405 //time_before_setVal=clock();
406 std::map< std::string, std::vector< BarlowCache > >::iterator iter_cache;
407 for( iter_cache = _barlowCache.begin(); iter_cache != _barlowCache.end(); ++iter_cache ) {
408
409 std::string channel_name = (*iter_cache).first;
410 std::vector< BarlowCache >& channel_cache = (*iter_cache).second;
411
412 /* Slower way to find the channel vector:
413 // Get the vector of bin uncertainty caches for this channel
414 if( _barlowCache.find( channel_name ) == _barlowCache.end() ) {
415 std::cout << "Error: channel: " << channel_name
416 << " not found in barlow Cache" << std::endl;
417 throw runtime_error("Channel not in barlow cache");
418 }
419
420 std::vector< BarlowCache >& channel_cache = _barlowCache[ channel_name ];
421 */
422
423 // Loop over the bins in the cache
424 // Set all gamma's to 0
425 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
426 BarlowCache& bin_cache = channel_cache.at(i);
427 if( !bin_cache.hasStatUncert ) continue;
428 RooRealVar* gamma = bin_cache.gamma;
429 gamma->setVal(0.0);
430 }
431 std::vector< double > nu_b_vec( channel_cache.size() );
432 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
433 BarlowCache& bin_cache = channel_cache.at(i);
434 if( !bin_cache.hasStatUncert ) continue;
435
436 RooAbsPdf* sum_pdf = static_cast<RooAbsPdf*>(bin_cache.sumPdf);
437 RooArgSet* obsSet = bin_cache.observables;
438 double binVolume = bin_cache.binVolume;
439
440 bin_cache.SetBinCenter();
441 double nu_b = sum_pdf->getVal(*obsSet)*sum_pdf->expectedEvents(*obsSet)*binVolume;
442 nu_b_vec.at(i) = nu_b;
443 }
444
445 // Loop over the bins in the cache
446 // Set all gamma's to 1
447 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
448 BarlowCache& bin_cache = channel_cache.at(i);
449 if( !bin_cache.hasStatUncert ) continue;
450 RooRealVar* gamma = bin_cache.gamma;
451 gamma->setVal(1.0);
452 }
453 std::vector< double > nu_b_stat_vec( channel_cache.size() );
454 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
455 BarlowCache& bin_cache = channel_cache.at(i);
456 if( !bin_cache.hasStatUncert ) continue;
457
458 RooAbsPdf* sum_pdf = static_cast<RooAbsPdf*>(bin_cache.sumPdf);
459 RooArgSet* obsSet = bin_cache.observables;
460 double binVolume = bin_cache.binVolume;
461
462 bin_cache.SetBinCenter();
463 double nu_b_stat = sum_pdf->getVal(*obsSet)*sum_pdf->expectedEvents(*obsSet)*binVolume - nu_b_vec.at(i);
464 nu_b_stat_vec.at(i) = nu_b_stat;
465 }
466 //time_after_setVal=clock();
467
468 // Done with the first loops.
469 // Now evaluating the function
470
471 //clock_t time_before_eval, time_after_eval;
472
473 // Loop over the bins in the cache
474 //time_before_eval=clock();
475 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
476
477 BarlowCache& bin_cache = channel_cache.at(i);
478
479 if( !bin_cache.hasStatUncert ) {
480 //std::cout << "Bin: " << i << " of " << channel_cache.size()
481 // << " in channel: " << channel_name
482 // << " doesn't have stat uncertainties" << std::endl;
483 continue;
484 }
485
486 // Set the observable to the bin center
487 bin_cache.SetBinCenter();
488
489 // Get the cached objects
490 RooRealVar* gamma = bin_cache.gamma;
491 RooRealVar* tau = bin_cache.tau;
492 RooAbsReal* pois_mean = bin_cache.nom_pois_mean;
493 //RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
494 //RooArgSet* obsSet = bin_cache.observables;
495 //double binVolume = bin_cache.binVolume;
496
497 // Get the values necessary for
498 // the analytic minimization
499 double nu_b = nu_b_vec.at(i);
500 double nu_b_stat = nu_b_stat_vec.at(i);
501
502 double tau_val = tau->getVal();
503 double nData = bin_cache.nData;
504 double m_val = pois_mean->getVal();
505
506 // Initialize the minimized value of gamma
507 double gamma_hat_hat = 1.0;
508
509 // Check that the quadratic term is > 0
510 if(nu_b_stat > 0.00000001) {
511
512 double A = nu_b_stat*nu_b_stat + tau_val*nu_b_stat;
513 double B = nu_b*tau_val + nu_b*nu_b_stat - nData*nu_b_stat - m_val*nu_b_stat;
514 double C = -1*m_val*nu_b;
515
516 double discrim = B*B-4*A*C;
517
518 if( discrim < 0 ) {
519 std::cout << "Warning: Discriminant (B*B - 4AC) < 0" << std::endl;
520 std::cout << "Warning: Taking B*B - 4*A*C == 0" << std::endl;
521 discrim=0;
522 //throw runtime_error("BarlowBeestonLL::evaluate() : B*B - 4AC < 0");
523 }
524 if( A <= 0 ) {
525 std::cout << "Warning: A <= 0" << std::endl;
526 throw runtime_error("BarlowBeestonLL::evaluate() : A < 0");
527 }
528
529 gamma_hat_hat = ( -1*B + std::sqrt(discrim) ) / (2*A);
530 }
531
532 // If the quadratic term is 0, we simply
533 // use a linear equation
534 else {
535 gamma_hat_hat = m_val/tau_val;
536 }
537
538 // Check for NAN
539 if( TMath::IsNaN(gamma_hat_hat) ) {
540 std::cout << "ERROR: gamma hat hat is NAN" << std::endl;
541 throw runtime_error("BarlowBeestonLL::evaluate() : gamma hat hat is NAN");
542 }
543
544 if( gamma_hat_hat <= 0 ) {
545 std::cout << "WARNING: gamma hat hat <= 0. Setting to 0" << std::endl;
546 gamma_hat_hat = 0;
547 }
548
549 /*
550 std::cout << "n: " << bin_cache.nData << " "
551 << "nu_stat: " << nu_b_stat << " "
552 << "nu: " << nu_b << " "
553 << "tau: " << tau->getVal() << " "
554 << "m: " << pois_mean->getVal() << " "
555 << "A: " << A << " "
556 << "B: " << B << " "
557 << "C: " << C << " "
558 << "gamma hat hat: " << gamma_hat_hat
559 << std::endl;
560 */
561
562 gamma->setVal( gamma_hat_hat );
563
564 }
565
566 //time_after_eval=clock();
567
568 //float time_setVal = ((float) time_after_setVal - (float) time_before_setVal) / ((float) CLOCKS_PER_SEC);
569 //float time_eval = ((float) time_after_eval - (float) time_before_eval) / ((float) CLOCKS_PER_SEC);
570
571 /*
572 std::cout << "Barlow timing for channel: " << channel_name
573 << " SetVal: " << time_setVal
574 << " Eval: " << time_eval
575 << std::endl;
576 */
577 }
578
579
580 return _nll;
581
582}
583
584
585
586/*
587////////////////////////////////////////////////////////////////////////////////
588/// Check that parameters and likelihood value for 'best fit' are still valid. If not,
589/// because the best fit has never been calculated, or because constant parameters have
590/// changed value or parameters have changed const/float status, the minimum is recalculated
591
592void RooStats::HistFactory::RooBarlowBeestonLL::validateAbsMin() const
593{
594 // Check if constant status of any of the parameters have changed
595 if (_absMinValid) {
596 _piter->Reset() ;
597 RooAbsArg* par ;
598 while((par=(RooAbsArg*)_piter->Next())) {
599 if (_paramFixed[par->GetName()] != par->isConstant()) {
600 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") constant status of parameter " << par->GetName() << " has changed from "
601 << (_paramFixed[par->GetName()]?"fixed":"floating") << " to " << (par->isConstant()?"fixed":"floating")
602 << ", recalculating absolute minimum" << endl ;
603 _absMinValid = false ;
604 break ;
605 }
606 }
607 }
608
609
610 // If we don't have the absolute minimum w.r.t all observables, calculate that first
611 if (!_absMinValid) {
612
613 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") determining minimum likelihood for current configurations w.r.t all observable" << endl ;
614
615
616 // Save current values of non-marginalized parameters
617 std::unique_ptr<RooArgSet> obsStart{(RooArgSet*) _obs.snapshot(false)};
618
619 // Start from previous global minimum
620 if (_paramAbsMin.size()>0) {
621 const_cast<RooSetProxy&>(_par).assignValueOnly(_paramAbsMin) ;
622 }
623 if (_obsAbsMin.size()>0) {
624 const_cast<RooSetProxy&>(_obs).assignValueOnly(_obsAbsMin) ;
625 }
626
627 // Find minimum with all observables floating
628 const_cast<RooSetProxy&>(_obs).setAttribAll("Constant",false) ;
629 _minuit->migrad() ;
630
631 // Save value and remember
632 _absMin = _nll ;
633 _absMinValid = true ;
634
635 // Save parameter values at abs minimum as well
636 _paramAbsMin.removeAll() ;
637
638 // Only store non-constant parameters here!
639 std::unique_ptr<RooArgSet> tmp{(RooArgSet*) _par.selectByAttrib("Constant",false)};
640 _paramAbsMin.addClone(*tmp) ;
641
642 _obsAbsMin.addClone(_obs) ;
643
644 // Save constant status of all parameters
645 _piter->Reset() ;
646 RooAbsArg* par ;
647 while((par=(RooAbsArg*)_piter->Next())) {
648 _paramFixed[par->GetName()] = par->isConstant() ;
649 }
650
651 if (dologI(Minimization)) {
652 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") minimum found at (" ;
653
654 RooAbsReal* arg ;
655 bool first=true ;
656 _oiter->Reset() ;
657 while ((arg=(RooAbsReal*)_oiter->Next())) {
658 ccxcoutI(Minimization) << (first?"":", ") << arg->GetName() << "=" << arg->getVal() ;
659 first=false ;
660 }
661 ccxcoutI(Minimization) << ")" << endl ;
662 }
663
664 // Restore original parameter values
665 const_cast<RooSetProxy&>(_obs) = *obsStart ;
666
667 }
668}
669*/
#define ClassImp(name)
Definition Rtypes.h:382
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t target
char name[80]
Definition TGX11.cxx:110
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
const RooArgSet * get(Int_t masterIdx) const
Int_t numBins() const
RooAbsReal & getParameter() const
double binVolume() const
bool isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:304
RooFit::OwningPtr< RooArgSet > getParameters(const RooAbsData *data, bool stripDisconnected=true) const
Create a list of leaf nodes in the arg tree starting with ourself as top node that don't match any of...
virtual bool remove(const RooAbsArg &var, bool silent=false, bool matchByNameOnly=false)
Remove the specified argument from our list.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
void reserve(Storage_t::size_type count)
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual double expectedEvents(const RooArgSet *nset) const
Return expected number of events to be used in calculation of extended likelihood.
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
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
Variable that can be changed from the outside.
Definition RooRealVar.h:37
Facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
Class RooBarlowBeestonLL implements the profile likelihood estimator for a given likelihood and set o...
double evaluate() const override
Optimized implementation of createProfile for profile likelihoods.
bool getParameters(const RooArgSet *depList, RooArgSet &outputSet, bool stripDisconnected=true) const override
Fills a list with leaf nodes in the arg tree starting with ourself as top node that don't match any o...
RooBarlowBeestonLL()
Default constructor. Should only be used by proof.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
bool getStatUncertaintyFromChannel(RooAbsPdf *channel, ParamHistFunc *&paramfunc, RooArgList *gammaList)
void FactorizeHistFactoryPdf(const RooArgSet &, RooAbsPdf &, RooArgList &, RooArgList &)
void getDataValuesForObservables(std::map< std::string, std::vector< double > > &ChannelBinDataMap, RooAbsData *data, RooAbsPdf *simPdf)
RooAbsPdf * getSumPdfFromChannel(RooAbsPdf *channel)
int getStatUncertaintyConstraintTerm(RooArgList *constraints, RooRealVar *gamma_stat, RooAbsReal *&pois_mean, RooRealVar *&tau)
Bool_t IsNaN(Double_t x)
Definition TMath.h:896