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
27#include "RooFit.h"
29#include "RooAbsReal.h"
30#include "RooAbsData.h"
31//#include "RooMinuit.h"
32#include "RooMsgService.h"
33#include "RooRealVar.h"
34#include "RooNLLVar.h"
35
37#include "RooProdPdf.h"
38#include "RooCategory.h"
39#include "RooSimultaneous.h"
40#include "RooArgList.h"
42
45
46using namespace std ;
47
49
50
51////////////////////////////////////////////////////////////////////////////////
52
54 RooAbsReal("RooBarlowBeestonLL","RooBarlowBeestonLL"),
55 _nll(),
56// _obs("paramOfInterest","Parameters of interest",this),
57// _par("nuisanceParam","Nuisance parameters",this,kFALSE,kFALSE),
58 _pdf(NULL), _data(NULL)
59{
60 // Default constructor
61 // Should only be used by proof.
62 // _piter = _par.createIterator() ;
63 // _oiter = _obs.createIterator() ;
64}
65
66
67////////////////////////////////////////////////////////////////////////////////
68
70 RooAbsReal& nllIn /*, const RooArgSet& observables*/) :
71 RooAbsReal(name,title),
72 _nll("input","-log(L) function",this,nllIn),
73 // _obs("paramOfInterest","Parameters of interest",this),
74 // _par("nuisanceParam","Nuisance parameters",this,kFALSE,kFALSE),
75 _pdf(NULL), _data(NULL)
76{
77 // Constructor of profile likelihood given input likelihood nll w.r.t
78 // the given set of variables. The input log likelihood is minimized w.r.t
79 // to all other variables of the likelihood at each evaluation and the
80 // value of the global log likelihood minimum is always subtracted.
81
82 // Determine actual parameters and observables
83 /*
84 RooArgSet* actualObs = nllIn.getObservables(observables) ;
85 RooArgSet* actualPars = nllIn.getParameters(observables) ;
86
87 _obs.add(*actualObs) ;
88 _par.add(*actualPars) ;
89
90 delete actualObs ;
91 delete actualPars ;
92
93 _piter = _par.createIterator() ;
94 _oiter = _obs.createIterator() ;
95 */
96}
97
98
99
100////////////////////////////////////////////////////////////////////////////////
101
103 RooAbsReal(other,name),
104 _nll("nll",this,other._nll),
105 // _obs("obs",this,other._obs),
106 // _par("par",this,other._par),
107 _pdf(NULL), _data(NULL),
108 _paramFixed(other._paramFixed)
109{
110 // Copy constructor
111
112 // _piter = _par.createIterator() ;
113 // _oiter = _obs.createIterator() ;
114
115 // _paramAbsMin.addClone(other._paramAbsMin) ;
116 // _obsAbsMin.addClone(other._obsAbsMin) ;
117
118}
119
120
121
122////////////////////////////////////////////////////////////////////////////////
123/// Destructor
124
126{
127 // Delete instance of minuit if it was ever instantiated
128 // if (_minuit) {
129 // delete _minuit ;
130 // }
131
132
133 //delete _piter ;
134 //delete _oiter ;
135}
136
137
138////////////////////////////////////////////////////////////////////////////////
139
141 TIterator* iter = bin_center->createIterator() ;
142 RooRealVar* var;
143 while((var=(RooRealVar*)iter->Next())) {
144 RooRealVar* target = (RooRealVar*) observables->find(var->GetName()) ;
145 target->setVal(var->getVal()) ;
146 }
147 delete iter;
148}
149
150
151////////////////////////////////////////////////////////////////////////////////
152
154 bool verbose=false;
155
156 if(!_data) {
157 std::cout << "Error: Must initialize data before initializing cache" << std::endl;
158 throw std::runtime_error("Uninitialized Data");
159 }
160 if(!_pdf) {
161 std::cout << "Error: Must initialize model pdf before initializing cache" << std::endl;
162 throw std::runtime_error("Uninitialized model pdf");
163 }
164
165 // Get the data bin values for all channels and bins
166 std::map< std::string, std::vector<double> > ChannelBinDataMap;
167 getDataValuesForObservables( ChannelBinDataMap, _data, _pdf );
168
169 // Get a list of constraint terms
170 RooArgList obsTerms;
171 RooArgList constraints;
172 RooArgSet* obsSet = _pdf->getObservables(*_data);
173 FactorizeHistFactoryPdf(*obsSet, *_pdf, obsTerms, constraints);
174
175 if( obsTerms.getSize() == 0 ) {
176 std::cout << "Error: Found no observable terms with pdf: " << _pdf->GetName()
177 << " using dataset: " << _data->GetName() << std::endl;
178 return;
179 }
180 if( constraints.getSize() == 0 ) {
181 std::cout << "Error: Found no constraint terms with pdf: " << _pdf->GetName()
182 << " using dataset: " << _data->GetName() << std::endl;
183 return;
184 }
185
186 /*
187 // Get the channels for this pdf
188 RooArgSet* channels = new RooArgSet();
189 RooArgSet* channelsWithConstraints = new RooArgSet();
190 getChannelsFromModel( _pdf, channels, channelsWithConstraints );
191 */
192
193 // Loop over the channels
194 RooSimultaneous* simPdf = (RooSimultaneous*) _pdf;
195 RooCategory* channelCat = (RooCategory*) (&simPdf->indexCat());
196 for (const auto& nameIdx : *channelCat) {
197
198 // Warning: channel cat name is not necesarily the same name
199 // as the pdf's (for example, if someone does edits)
200 RooAbsPdf* channelPdf = simPdf->getPdf(nameIdx.first.c_str());
201 std::string channel_name = channelPdf->GetName();
202
203 // First, we check if this channel uses Stat Uncertainties:
204 RooArgList* gammas = new RooArgList();
205 ParamHistFunc* param_func=NULL;
206 bool hasStatUncert = getStatUncertaintyFromChannel( channelPdf, param_func, gammas );
207 if( ! hasStatUncert ) {
208 if(verbose) {
209 std::cout << "Channel: " << channel_name
210 << " doesn't have statistical uncertainties"
211 << std::endl;
212 }
213 continue;
214 }
215 else {
216 if(verbose) std::cout << "Found ParamHistFunc: " << param_func->GetName() << std::endl;
217 }
218
219 // Now, loop over the bins in this channel
220 // To Do: Check that the index convention
221 // still works for 2-d (ie matches the
222 // convention in ParamHistFunc, etc)
223 int num_bins = param_func->numBins();
224
225 // Initialize the vector to the number of bins, allowing
226 // us to skip gamma's if they are constant
227
228 std::vector<BarlowCache> temp_cache( num_bins );
229 bool channel_has_stat_uncertainty=false;
230
231 for( Int_t bin_index = 0; bin_index < num_bins; ++bin_index ) {
232
233 // Create a cache object
234 BarlowCache cache;
235
236 // Get the gamma for this bin, and skip if it's constant
237 RooRealVar* gamma_stat = &(param_func->getParameter(bin_index));
238 if( gamma_stat->isConstant() ) {
239 if(verbose) std::cout << "Ignoring constant gamma: " << gamma_stat->GetName() << std::endl;
240 continue;
241 }
242 else {
243 cache.hasStatUncert=true;
244 channel_has_stat_uncertainty=true;
245 cache.gamma = gamma_stat;
246 _statUncertParams.insert( gamma_stat->GetName() );
247 }
248
249 // Store a snapshot of the bin center
250 RooArgSet* bin_center = (RooArgSet*) param_func->get( bin_index )->snapshot();
251 cache.bin_center = bin_center;
252 cache.observables = obsSet;
253
254 cache.binVolume = param_func->binVolume();
255
256 // Delete this part, simply a comment
257 RooArgList obs_list( *(cache.bin_center) );
258
259 // Get the gamma's constraint term
260 RooAbsReal* pois_mean = NULL;
261 RooRealVar* tau = NULL;
262 getStatUncertaintyConstraintTerm( &constraints, gamma_stat, pois_mean, tau );
263 if( !tau || !pois_mean ) {
264 std::cout << "Failed to find pois mean or tau parameter for " << gamma_stat->GetName() << std::endl;
265 }
266 else {
267 if(verbose) std::cout << "Found pois mean and tau for parameter: " << gamma_stat->GetName()
268 << " tau: " << tau->GetName() << " " << tau->getVal()
269 << " pois_mean: " << pois_mean->GetName() << " " << pois_mean->getVal()
270 << std::endl;
271 }
272
273 cache.tau = tau;
274 cache.nom_pois_mean = pois_mean;
275
276 // Get the RooRealSumPdf
277 RooAbsPdf* sum_pdf = getSumPdfFromChannel( channelPdf );
278 if( sum_pdf == NULL ) {
279 std::cout << "Failed to find RooRealSumPdf in channel " << channel_name
280 << ", therefor skipping this channel for analytic uncertainty minimization"
281 << std::endl;
282 channel_has_stat_uncertainty=false;
283 break;
284 }
285 cache.sumPdf = sum_pdf;
286
287 // And set the data value for this bin
288 if( ChannelBinDataMap.find(channel_name) == ChannelBinDataMap.end() ) {
289 std::cout << "Error: channel with name: " << channel_name
290 << " not found in BinDataMap" << std::endl;
291 throw runtime_error("BinDataMap");
292 }
293 double nData = ChannelBinDataMap[channel_name].at(bin_index);
294 cache.nData = nData;
295
296 temp_cache.at(bin_index) = cache;
297 //_barlowCache[channel_name].at(bin_index) = cache;
298
299 } // End: Loop over bins
300
301 if( channel_has_stat_uncertainty ) {
302 std::cout << "Adding channel: " << channel_name
303 << " to the barlow cache" << std::endl;
304 _barlowCache[channel_name] = temp_cache;
305 }
306
307
308 } // End: Loop over channels
309
310
311
312 // Successfully initialized the cache
313 // Printing some info
314 /*
315 std::map< std::string, std::vector< BarlowCache > >::iterator iter_cache;
316 for( iter_cache = _barlowCache.begin(); iter_cache != _barlowCache.end(); ++iter_cache ) {
317
318 std::string channel_name = (*iter_cache).first;
319 std::vector< BarlowCache >& channel_cache = (*iter_cache).second;
320
321
322 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
323 BarlowCache& bin_cache = channel_cache.at(i);
324
325
326 RooRealVar* gamma = bin_cache.gamma;
327 RooRealVar* tau = bin_cache.tau;
328 RooAbsReal* pois_mean = bin_cache.nom_pois_mean;
329 RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
330 double binVolume = bin_cache.binVolume;
331
332
333 if( !bin_cache.hasStatUncert ) {
334 std::cout << "Barlow Cache for Channel: " << channel_name
335 << " Bin: " << i
336 << " NO STAT UNCERT"
337 << std::endl;
338 }
339 else {
340 std::cout << "Barlow Cache for Channel: " << channel_name
341 << " Bin: " << i
342 << " gamma: " << gamma->GetName()
343 << " tau: " << tau->GetName()
344 << " pois_mean: " << pois_mean->GetName()
345 << " sum_pdf: " << sum_pdf->GetName()
346 << " binVolume: " << binVolume
347 << std::endl;
348 }
349
350 }
351 }
352 */
353
354}
355
356
357////////////////////////////////////////////////////////////////////////////////
358
360 RooArgSet& outputSet,
361 bool stripDisconnected) const {
362 bool errorInBaseCall = RooAbsArg::getParameters( depList, outputSet, stripDisconnected );
363
364 for (auto const& arg : outputSet) {
365
366 // If there is a gamma in the name,
367 // strip it from the list of dependencies
368
369 if( _statUncertParams.find(arg->GetName()) != _statUncertParams.end() ) {
370 outputSet.remove( *arg, true );
371 }
372
373 }
374
375 return errorInBaseCall || false;
376
377}
378
379
380/*
381////////////////////////////////////////////////////////////////////////////////
382
383const RooArgSet& RooStats::HistFactory::RooBarlowBeestonLL::bestFitParams() const
384{
385 validateAbsMin() ;
386 return _paramAbsMin ;
387}
388
389
390////////////////////////////////////////////////////////////////////////////////
391
392const RooArgSet& RooStats::HistFactory::RooBarlowBeestonLL::bestFitObs() const
393{
394 validateAbsMin() ;
395 return _obsAbsMin ;
396}
397*/
398
399
400
401////////////////////////////////////////////////////////////////////////////////
402/// Optimized implementation of createProfile for profile likelihoods.
403/// Return profile of original function in terms of stated parameters
404/// of interest rather than profiling recursively.
405
406/*
407RooAbsReal* RooStats::HistFactory::RooBarlowBeestonLL::createProfile(const RooArgSet& paramsOfInterest)
408{
409 return nll().createProfile(paramsOfInterest) ;
410}
411*/
412
413
414/*
415void RooStats::HistFactory::RooBarlowBeestonLL::FactorizePdf(const RooArgSet &observables, RooAbsPdf &pdf, RooArgList &obsTerms, RooArgList &constraints) const {
416 // utility function to factorize constraint terms from a pdf
417 // (from G. Petrucciani)
418 const std::type_info & id = typeid(pdf);
419 if (id == typeid(RooProdPdf)) {
420 RooProdPdf *prod = dynamic_cast<RooProdPdf *>(&pdf);
421 RooArgList list(prod->pdfList());
422 for (int i = 0, n = list.getSize(); i < n; ++i) {
423 RooAbsPdf *pdfi = (RooAbsPdf *) list.at(i);
424 FactorizePdf(observables, *pdfi, obsTerms, constraints);
425 }
426 } else if (id == typeid(RooSimultaneous) ) { //|| id == typeid(RooSimultaneousOpt)) {
427 RooSimultaneous *sim = dynamic_cast<RooSimultaneous *>(&pdf);
428 RooAbsCategoryLValue *cat = (RooAbsCategoryLValue *) sim->indexCat().Clone();
429 for (int ic = 0, nc = cat->numBins((const char *)0); ic < nc; ++ic) {
430 cat->setBin(ic);
431 FactorizePdf(observables, *sim->getPdf(cat->getLabel()), obsTerms, constraints);
432 }
433 delete cat;
434 } else if (pdf.dependsOn(observables)) {
435 if (!obsTerms.contains(pdf)) obsTerms.add(pdf);
436 } else {
437 if (!constraints.contains(pdf)) constraints.add(pdf);
438 }
439}
440*/
441
442
443
444////////////////////////////////////////////////////////////////////////////////
445
447{
448 /*
449 // Loop over the cached bins and channels
450 RooArgSet* channels = new RooArgSet();
451 RooArgSet* channelsWithConstraints = new RooArgSet();
452 RooStats::getChannelsFromModel( _pdf, channels, channelsWithConstraints );
453
454 // Loop over channels
455 TIterator* iter_channels = channelsWithConstraints->createIterator();
456 RooAbsPdf* channelPdf=NULL;
457 while(( channelPdf=(RooAbsPdf*)iter_channels->Next() )) {
458 std::string channel_name = channelPdf->GetName(); //RooStats::channelNameFromPdf( channelPdf );
459 */
460
461 // Loop over the channels (keys to the map)
462 //clock_t time_before_setVal, time_after_setVal;
463 //time_before_setVal=clock();
464 std::map< std::string, std::vector< BarlowCache > >::iterator iter_cache;
465 for( iter_cache = _barlowCache.begin(); iter_cache != _barlowCache.end(); ++iter_cache ) {
466
467 std::string channel_name = (*iter_cache).first;
468 std::vector< BarlowCache >& channel_cache = (*iter_cache).second;
469
470 /* Slower way to find the channel vector:
471 // Get the vector of bin uncertainty caches for this channel
472 if( _barlowCache.find( channel_name ) == _barlowCache.end() ) {
473 std::cout << "Error: channel: " << channel_name
474 << " not found in barlow Cache" << std::endl;
475 throw runtime_error("Channel not in barlow cache");
476 }
477
478 std::vector< BarlowCache >& channel_cache = _barlowCache[ channel_name ];
479 */
480
481 // Loop over the bins in the cache
482 // Set all gamma's to 0
483 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
484 BarlowCache& bin_cache = channel_cache.at(i);
485 if( !bin_cache.hasStatUncert ) continue;
486 RooRealVar* gamma = bin_cache.gamma;
487 gamma->setVal(0.0);
488 }
489 std::vector< double > nu_b_vec( channel_cache.size() );
490 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
491 BarlowCache& bin_cache = channel_cache.at(i);
492 if( !bin_cache.hasStatUncert ) continue;
493
494 RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
495 RooArgSet* obsSet = bin_cache.observables;
496 double binVolume = bin_cache.binVolume;
497
498 bin_cache.SetBinCenter();
499 double nu_b = sum_pdf->getVal(*obsSet)*sum_pdf->expectedEvents(*obsSet)*binVolume;
500 nu_b_vec.at(i) = nu_b;
501 }
502
503 // Loop over the bins in the cache
504 // Set all gamma's to 1
505 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
506 BarlowCache& bin_cache = channel_cache.at(i);
507 if( !bin_cache.hasStatUncert ) continue;
508 RooRealVar* gamma = bin_cache.gamma;
509 gamma->setVal(1.0);
510 }
511 std::vector< double > nu_b_stat_vec( channel_cache.size() );
512 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
513 BarlowCache& bin_cache = channel_cache.at(i);
514 if( !bin_cache.hasStatUncert ) continue;
515
516 RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
517 RooArgSet* obsSet = bin_cache.observables;
518 double binVolume = bin_cache.binVolume;
519
520 bin_cache.SetBinCenter();
521 double nu_b_stat = sum_pdf->getVal(*obsSet)*sum_pdf->expectedEvents(*obsSet)*binVolume - nu_b_vec.at(i);
522 nu_b_stat_vec.at(i) = nu_b_stat;
523 }
524 //time_after_setVal=clock();
525
526 // Done with the first loops.
527 // Now evaluating the function
528
529 //clock_t time_before_eval, time_after_eval;
530
531 // Loop over the bins in the cache
532 //time_before_eval=clock();
533 for( unsigned int i = 0; i < channel_cache.size(); ++i ) {
534
535 BarlowCache& bin_cache = channel_cache.at(i);
536
537 if( !bin_cache.hasStatUncert ) {
538 //std::cout << "Bin: " << i << " of " << channel_cache.size()
539 // << " in channel: " << channel_name
540 // << " doesn't have stat uncertainties" << std::endl;
541 continue;
542 }
543
544 // Set the observable to the bin center
545 bin_cache.SetBinCenter();
546
547 // Get the cached objects
548 RooRealVar* gamma = bin_cache.gamma;
549 RooRealVar* tau = bin_cache.tau;
550 RooAbsReal* pois_mean = bin_cache.nom_pois_mean;
551 //RooAbsPdf* sum_pdf = (RooAbsPdf*) bin_cache.sumPdf;
552 //RooArgSet* obsSet = bin_cache.observables;
553 //double binVolume = bin_cache.binVolume;
554
555 // Get the values necessary for
556 // the analytic minimization
557 double nu_b = nu_b_vec.at(i);
558 double nu_b_stat = nu_b_stat_vec.at(i);
559
560 double tau_val = tau->getVal();
561 double nData = bin_cache.nData;
562 double m_val = pois_mean->getVal();
563
564 // Initialize the minimized value of gamma
565 double gamma_hat_hat = 1.0;
566
567 // Check that the quadratic term is > 0
568 if(nu_b_stat > 0.00000001) {
569
570 double A = nu_b_stat*nu_b_stat + tau_val*nu_b_stat;
571 double B = nu_b*tau_val + nu_b*nu_b_stat - nData*nu_b_stat - m_val*nu_b_stat;
572 double C = -1*m_val*nu_b;
573
574 double discrim = B*B-4*A*C;
575
576 if( discrim < 0 ) {
577 std::cout << "Warning: Discriminant (B*B - 4AC) < 0" << std::endl;
578 std::cout << "Warning: Taking B*B - 4*A*C == 0" << std::endl;
579 discrim=0;
580 //throw runtime_error("BarlowBeestonLL::evaluate() : B*B - 4AC < 0");
581 }
582 if( A <= 0 ) {
583 std::cout << "Warning: A <= 0" << std::endl;
584 throw runtime_error("BarlowBeestonLL::evaluate() : A < 0");
585 }
586
587 gamma_hat_hat = ( -1*B + TMath::Sqrt(discrim) ) / (2*A);
588 }
589
590 // If the quadratic term is 0, we simply
591 // use a linear equation
592 else {
593 gamma_hat_hat = m_val/tau_val;
594 }
595
596 // Check for NAN
597 if( TMath::IsNaN(gamma_hat_hat) ) {
598 std::cout << "ERROR: gamma hat hat is NAN" << std::endl;
599 throw runtime_error("BarlowBeestonLL::evaluate() : gamma hat hat is NAN");
600 }
601
602 if( gamma_hat_hat <= 0 ) {
603 std::cout << "WARNING: gamma hat hat <= 0. Setting to 0" << std::endl;
604 gamma_hat_hat = 0;
605 }
606
607 /*
608 std::cout << "n: " << bin_cache.nData << " "
609 << "nu_stat: " << nu_b_stat << " "
610 << "nu: " << nu_b << " "
611 << "tau: " << tau->getVal() << " "
612 << "m: " << pois_mean->getVal() << " "
613 << "A: " << A << " "
614 << "B: " << B << " "
615 << "C: " << C << " "
616 << "gamma hat hat: " << gamma_hat_hat
617 << std::endl;
618 */
619
620 gamma->setVal( gamma_hat_hat );
621
622 }
623
624 //time_after_eval=clock();
625
626 //float time_setVal = ((float) time_after_setVal - (float) time_before_setVal) / ((float) CLOCKS_PER_SEC);
627 //float time_eval = ((float) time_after_eval - (float) time_before_eval) / ((float) CLOCKS_PER_SEC);
628
629 /*
630 std::cout << "Barlow timing for channel: " << channel_name
631 << " SetVal: " << time_setVal
632 << " Eval: " << time_eval
633 << std::endl;
634 */
635 }
636
637
638 return _nll;
639
640}
641
642
643
644/*
645////////////////////////////////////////////////////////////////////////////////
646/// Check that parameters and likelihood value for 'best fit' are still valid. If not,
647/// because the best fit has never been calculated, or because constant parameters have
648/// changed value or parameters have changed const/float status, the minimum is recalculated
649
650void RooStats::HistFactory::RooBarlowBeestonLL::validateAbsMin() const
651{
652 // Check if constant status of any of the parameters have changed
653 if (_absMinValid) {
654 _piter->Reset() ;
655 RooAbsArg* par ;
656 while((par=(RooAbsArg*)_piter->Next())) {
657 if (_paramFixed[par->GetName()] != par->isConstant()) {
658 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") constant status of parameter " << par->GetName() << " has changed from "
659 << (_paramFixed[par->GetName()]?"fixed":"floating") << " to " << (par->isConstant()?"fixed":"floating")
660 << ", recalculating absolute minimum" << endl ;
661 _absMinValid = kFALSE ;
662 break ;
663 }
664 }
665 }
666
667
668 // If we don't have the absolute minimum w.r.t all observables, calculate that first
669 if (!_absMinValid) {
670
671 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") determining minimum likelihood for current configurations w.r.t all observable" << endl ;
672
673
674 // Save current values of non-marginalized parameters
675 RooArgSet* obsStart = (RooArgSet*) _obs.snapshot(kFALSE) ;
676
677 // Start from previous global minimum
678 if (_paramAbsMin.getSize()>0) {
679 const_cast<RooSetProxy&>(_par).assignValueOnly(_paramAbsMin) ;
680 }
681 if (_obsAbsMin.getSize()>0) {
682 const_cast<RooSetProxy&>(_obs).assignValueOnly(_obsAbsMin) ;
683 }
684
685 // Find minimum with all observables floating
686 const_cast<RooSetProxy&>(_obs).setAttribAll("Constant",kFALSE) ;
687 _minuit->migrad() ;
688
689 // Save value and remember
690 _absMin = _nll ;
691 _absMinValid = kTRUE ;
692
693 // Save parameter values at abs minimum as well
694 _paramAbsMin.removeAll() ;
695
696 // Only store non-constant parameters here!
697 RooArgSet* tmp = (RooArgSet*) _par.selectByAttrib("Constant",kFALSE) ;
698 _paramAbsMin.addClone(*tmp) ;
699 delete tmp ;
700
701 _obsAbsMin.addClone(_obs) ;
702
703 // Save constant status of all parameters
704 _piter->Reset() ;
705 RooAbsArg* par ;
706 while((par=(RooAbsArg*)_piter->Next())) {
707 _paramFixed[par->GetName()] = par->isConstant() ;
708 }
709
710 if (dologI(Minimization)) {
711 cxcoutI(Minimization) << "RooStats::HistFactory::RooBarlowBeestonLL::evaluate(" << GetName() << ") minimum found at (" ;
712
713 RooAbsReal* arg ;
714 Bool_t first=kTRUE ;
715 _oiter->Reset() ;
716 while ((arg=(RooAbsReal*)_oiter->Next())) {
717 ccxcoutI(Minimization) << (first?"":", ") << arg->GetName() << "=" << arg->getVal() ;
718 first=kFALSE ;
719 }
720 ccxcoutI(Minimization) << ")" << endl ;
721 }
722
723 // Restore original parameter values
724 const_cast<RooSetProxy&>(_obs) = *obsStart ;
725 delete obsStart ;
726
727 }
728}
729*/
730
731
732////////////////////////////////////////////////////////////////////////////////
733
735 Bool_t /*nameChange*/, Bool_t /*isRecursive*/)
736{
737 /*
738 if (_minuit) {
739 delete _minuit ;
740 _minuit = 0 ;
741 }
742 */
743 return kFALSE ;
744}
745
746
const Bool_t kFALSE
Definition RtypesCore.h:92
#define ClassImp(name)
Definition Rtypes.h:364
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
RooRealVar & getParameter() const
Int_t numBins() const
double binVolume() const
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...
Bool_t isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:380
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
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 expectedEvents(const RooArgSet *nset) const
Return expected number of events from this p.d.f for use in extended likelihood calculations.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:21
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
RooArgSet * snapshot(bool deepCopy=true) const
Use RooAbsCollection::snapshot(), but return as RooArgSet.
Definition RooArgSet.h:118
RooCategory is an object to represent discrete states.
Definition RooCategory.h:27
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual void setVal(Double_t value)
Set value of variable to 'value'.
RooSimultaneous facilitates simultaneous fitting of multiple PDFs to subsets of a given dataset.
const RooAbsCategoryLValue & indexCat() const
RooAbsPdf * getPdf(const char *catName) const
Return the p.d.f associated with the given index category name.
Class RooBarlowBeestonLL implements the profile likelihood estimator for a given likelihood and set o...
bool getParameters(const RooArgSet *depList, RooArgSet &outputSet, bool stripDisconnected=true) const
Fills a list with leaf nodes in the arg tree starting with ourself as top node that don't match any o...
virtual Bool_t redirectServersHook(const RooAbsCollection &, Bool_t, Bool_t, Bool_t)
Function that is called at the end of redirectServers().
Double_t evaluate() const
Optimized implementation of createProfile for profile likelihoods.
Iterator abstract base class.
Definition TIterator.h:30
virtual TObject * Next()=0
virtual const char * GetName() const
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:892
Double_t Sqrt(Double_t x)
Definition TMath.h:691