CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BestFitSigmaTestStat.cc
Go to the documentation of this file.
1 #include "../interface/ProfiledLikelihoodRatioTestStatExt.h"
2 #include "../interface/BestFitSigmaTestStat.h"
3 #include "../interface/CascadeMinimizer.h"
4 #include "../interface/CloseCoutSentry.h"
5 #include "../interface/utils.h"
6 #include <stdexcept>
7 #include <RooRealVar.h>
8 #include "../interface/RooMinimizerOpt.h"
9 #include <RooFitResult.h>
10 #include <RooSimultaneous.h>
11 #include <RooCategory.h>
12 #include <RooRandom.h>
13 #include <Math/MinimizerOptions.h>
14 #include <RooStats/RooStatsUtils.h>
15 #include "../interface/ProfilingTools.h"
16 
17 //---- Uncomment this and run with --perfCounters to get statistics of successful and failed fits
18 #define DEBUG_FIT_STATUS
19 #ifdef DEBUG_FIT_STATUS
20 #define COUNT_ONE(x) PerfCounter::add(x);
21 #else
22 #define COUNT_ONE(X)
23 #endif
24 
25 //---- Uncomment this and set some of these to 1 to get more debugging
26 #if 0
27 #define DBG(X,Z) if (X) { Z; }
28 #define DBGV(X,Z) if (X>1) { Z; }
29 #define DBG_TestStat_params 0 // ProfiledLikelihoodRatioTestStatOpt::Evaluate; 1 = dump nlls; 2 = dump params at each eval
30 #define DBG_TestStat_NOFIT 0 // FIXME HACK: if set, don't profile the likelihood, just evaluate it
31 #define DBG_PLTestStat_ctor 0 // dump parameters in c-tor
32 #define DBG_PLTestStat_pars 0 // dump parameters in eval
33 #define DBG_PLTestStat_fit 0 // dump fit result
34 #define DBG_PLTestStat_main 1 // limited debugging (final NLLs, failed fits)
35 #else
36 #define DBG(X,Z)
37 #define DBGV(X,Z)
38 #endif
39 
40 
41 //============================================================================================================================================
43  const RooArgSet & observables,
44  RooAbsPdf &pdf,
45  const RooArgSet *nuisances,
46  const RooArgSet & params,
47  int verbosity)
48 :
49  pdf_(&pdf),
50  verbosity_(verbosity)
51 {
52  params.snapshot(snap_,false);
53  ((RooRealVar*)snap_.find(params.first()->GetName()))->setConstant(false);
54  poi_.add(snap_);
55  if (nuisances) { nuisances_.add(*nuisances); snap_.addClone(*nuisances, /*silent=*/true); }
56  params_.reset(pdf_->getParameters(observables));
57 }
58 
59 
60 Double_t BestFitSigmaTestStat::Evaluate(RooAbsData& data, RooArgSet& /*nullPOI*/)
61 {
62  // Take snapshot of initial state, to restore it at the end
63  RooArgSet initialState; params_->snapshot(initialState);
64 
65  // Initialize parameters
66  *params_ = snap_;
67 
68  // Initialize signal strength
69  RooRealVar *rIn = (RooRealVar *) poi_.first();
70  RooRealVar *r = (RooRealVar *) params_->find(rIn->GetName());
71  bool canKeepNLL = createNLL(*pdf_, data);
72 
73  double initialR = rIn->getVal();
74 
75  // Perform unconstrained minimization
76  r->setVal(initialR == 0 ? 0.5 : 0.5*initialR); //best guess
77  r->setConstant(false);
78 
79  std::cout << "Doing a fit for with " << r->GetName() << " in range [ " << r->getMin() << " , " << r->getMax() << "]" << std::endl;
80  std::cout << "Starting point is " << r->GetName() << " = " << r->getVal() << std::endl;
81  minNLL(/*constrained=*/false, r);
82  double bestFitR = r->getVal();
83  std::cout << "Fit result was " << r->GetName() << " = " << r->getVal() << std::endl;
84 
85  //Restore initial state, to avoid issues with ToyMCSampler
86  *params_ = initialState;
87 
88  if (!canKeepNLL) nll_.reset();
89 
90  return bestFitR;
91 }
92 
93 bool BestFitSigmaTestStat::createNLL(RooAbsPdf &pdf, RooAbsData &data)
94 {
95  if (typeid(pdf) == typeid(RooSimultaneousOpt)) {
96  if (nll_.get() == 0) nll_.reset(pdf.createNLL(data, RooFit::Constrain(nuisances_)));
97  else ((cacheutils::CachingSimNLL&)(*nll_)).setData(data);
98  return true;
99  } else {
100  nll_.reset(pdf.createNLL(data, RooFit::Constrain(nuisances_)));
101  return false;
102  }
103 }
104 
105 double BestFitSigmaTestStat::minNLL(bool constrained, RooRealVar *r)
106 {
108  CascadeMinimizer minim(*nll_, mode, r);
109  minim.minimize(verbosity_-2);
110  return nll_->getVal();
111 }
bool minimize(int verbose=0, bool cascade=true)
bool createNLL(RooAbsPdf &pdf, RooAbsData &data)
BestFitSigmaTestStat(const RooArgSet &observables, RooAbsPdf &pdf, const RooArgSet *nuisances, const RooArgSet &params, int verbosity=0)
std::auto_ptr< RooArgSet > params_
virtual Double_t Evaluate(RooAbsData &data, RooArgSet &nullPOI)
std::auto_ptr< RooAbsReal > nll_
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
tuple cout
Definition: gather_cfg.py:121
double minNLL(bool constrained, RooRealVar *r=0)