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