Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ProfileInspector.cxx
Go to the documentation of this file.
1// @(#)root/roostats:$Id: ProfileInspector.cxx 34109 2010-06-24 15:00:16Z moneta $
2
3/*************************************************************************
4 * Project: RooStats *
5 * Package: RooFit/RooStats *
6 * Authors: *
7 * Kyle Cranmer, Lorenzo Moneta, Gregory Schott, Wouter Verkerke *
8 * Akira Shibata
9 *************************************************************************
10 * Copyright (C) 1995-2008, Rene Brun and Fons Rademakers. *
11 * All rights reserved. *
12 * *
13 * For the licensing terms see $ROOTSYS/LICENSE. *
14 * For the list of contributors see $ROOTSYS/README/CREDITS. *
15 *************************************************************************/
16
17
18/** \class RooStats::ProfileInspector
19 \ingroup Roostats
20
21Utility class to plot conditional MLE of nuisance parameters vs. Parameters of Interest
22
23*/
24
25
27#include "RooRealVar.h"
28#include "RooAbsReal.h"
29#include "RooArgSet.h"
30#include "RooAbsPdf.h"
31#include "RooCurve.h"
32#include "TAxis.h"
33
34
35using namespace RooStats;
36using std::string, std::map;
37
38////////////////////////////////////////////////////////////////////////////////
39
43
44////////////////////////////////////////////////////////////////////////////////
45/// ProfileInspector destructor
46
50
51////////////////////////////////////////////////////////////////////////////////
52///
53/// This tool makes a plot of the conditional maximum likelihood estimate of the nuisance parameter
54/// vs the parameter of interest
55///
56/// This enables you to discover if any of the nuisance parameters are behaving strangely
57/// curve is the optional parameters, when used you can specify the points previously scanned
58/// in the process of plotOn or createHistogram.
59/// To do this, you can do the following after the plot has been made:
60/// ~~~ {.cpp}
61/// profile, RooRealVar * poi, RooCurve * curve ){
62///RooCurve * curve = 0;
63/// ~~~
64
66{
67
68 const RooArgSet* poi_set = config->GetParametersOfInterest();
70 RooAbsPdf* pdf = config->GetPdf();
71
72
73 if(!poi_set){
74 std::cout << "no parameters of interest" << std::endl;
75 return nullptr;
76 }
77
78 if(poi_set->size()!=1){
79 std::cout << "only one parameter of interest is supported currently" << std::endl;
80 return nullptr;
81 }
82 RooRealVar* poi = static_cast<RooRealVar*>(poi_set->first());
83
84
85 if(!nuis_params){
86 std::cout << "no nuisance parameters" << std::endl;
87 return nullptr;
88 }
89
90 if(!pdf){
91 std::cout << "pdf not set" << std::endl;
92 return nullptr;
93 }
94
95 std::unique_ptr<RooAbsReal> nll{pdf->createNLL(data)};
96 std::unique_ptr<RooAbsReal> profile{nll->createProfile(*poi)};
97
98 TList * list = new TList;
99 Int_t curve_N=100;
100 double* curve_x=nullptr;
101// if(curve){
102// curve_N=curve->GetN();
103// curve_x=curve->GetX();
104// } else {
105 double max = dynamic_cast<RooAbsRealLValue*>(poi)->getMax();
106 double min = dynamic_cast<RooAbsRealLValue*>(poi)->getMin();
107 double step = (max-min)/(curve_N-1);
108 curve_x=new double[curve_N];
109 for(int i=0; i<curve_N; ++i){
110 curve_x[i]=min+step*i;
111 }
112// }
113
115 for(int i=0; i<curve_N; i++){
116 poi->setVal(curve_x[i]);
117 profile->getVal();
118
120 if(nuis_param) {
121 string name = nuis_param->GetName();
122 if(nuis_params->empty()) continue;
123 if(nuis_param && (! nuis_param->isConstant())){
124 if(name_val.find(name)==name_val.end()) name_val[name]=std::vector<double>(curve_N);
125 name_val[name][i]=nuis_param->getVal();
126
127 if(i==curve_N-1){
128 TGraph* g = new TGraph(curve_N, curve_x, &(name_val[name].front()));
129 g->SetName((name+"_"+string(poi->GetName())+"_profile").c_str());
130 g->GetXaxis()->SetTitle(poi->GetName());
131 g->GetYaxis()->SetTitle(nuis_param->GetName());
132 g->SetTitle("");
133 list->Add(g);
134 }
135 }
136 }
137 }
138 }
139
140 delete [] curve_x;
141
142 return list;
143}
#define g(i)
Definition RSha256.hxx:105
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
const_iterator end() const
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
RooFit::OwningPtr< RooAbsReal > createNLL(RooAbsData &data, CmdArgs_t const &... cmdArgs)
Construct representation of -log(L) of PDF with given dataset.
Definition RooAbsPdf.h:163
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition ModelConfig.h:35
const RooArgSet * GetParametersOfInterest() const
get RooArgSet containing the parameter of interest (return nullptr if not existing)
const RooArgSet * GetNuisanceParameters() const
get RooArgSet containing the nuisance parameters (return nullptr if not existing)
RooAbsPdf * GetPdf() const
get model PDF (return nullptr if pdf has not been specified or does not exist)
TList * GetListOfProfilePlots(RooAbsData &data, RooStats::ModelConfig *config)
This tool makes a plot of the conditional maximum likelihood estimate of the nuisance parameter vs th...
virtual ~ProfileInspector()
ProfileInspector destructor.
A TGraph is an object made of two arrays X and Y with npoints each.
Definition TGraph.h:41
A doubly linked list.
Definition TList.h:38
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Namespace for the RooStats classes.
Definition CodegenImpl.h:58