Logo ROOT   6.30.04
Reference Guide
 All Namespaces Files Pages
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 
21 Utility 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 "RooArgSet.h"
32 #include "RooCurve.h"
33 #include "TAxis.h"
34 
35 ClassImp(RooStats::ProfileInspector);
36 
37 using namespace RooStats;
38 using namespace std;
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
42 ProfileInspector::ProfileInspector()
43 {
44 }
45 
46 ////////////////////////////////////////////////////////////////////////////////
47 /// ProfileInspector destructor
48 
49 ProfileInspector::~ProfileInspector()
50 {
51 }
52 
53 ////////////////////////////////////////////////////////////////////////////////
54 ///
55 /// This tool makes a plot of the conditional maximum likelihood estimate of the nuisance parameter
56 /// vs the parameter of interest
57 ///
58 /// This enables you to discover if any of the nuisance parameters are behaving strangely
59 /// curve is the optional parameters, when used you can specify the points previously scanned
60 /// in the process of plotOn or createHistogram.
61 /// To do this, you can do the following after the plot has been made:
62 /// ~~~ {.cpp}
63 /// profile, RooRealVar * poi, RooCurve * curve ){
64 ///RooCurve * curve = 0;
65 /// ~~~
66 
67 TList* ProfileInspector::GetListOfProfilePlots( RooAbsData& data, RooStats::ModelConfig * config)
68 {
69 
70  const RooArgSet* poi_set = config->GetParametersOfInterest();
71  const RooArgSet* nuis_params=config->GetNuisanceParameters();
72  RooAbsPdf* pdf = config->GetPdf();
73 
74 
75  if(!poi_set){
76  cout << "no parameters of interest" << endl;
77  return 0;
78  }
79 
80  if(poi_set->getSize()!=1){
81  cout << "only one parameter of interest is supported currently" << endl;
82  return 0;
83  }
84  RooRealVar* poi = (RooRealVar*) poi_set->first();
85 
86 
87  if(!nuis_params){
88  cout << "no nuisance parameters" << endl;
89  return 0;
90  }
91 
92  if(!pdf){
93  cout << "pdf not set" << endl;
94  return 0;
95  }
96 
97  RooAbsReal* nll = pdf->createNLL(data);
98  RooAbsReal* profile = nll->createProfile(*poi);
99 
100  TList * list = new TList;
101  Int_t curve_N=100;
102  Double_t* curve_x=0;
103 // if(curve){
104 // curve_N=curve->GetN();
105 // curve_x=curve->GetX();
106 // } else {
107  Double_t max = dynamic_cast<RooAbsRealLValue*>(poi)->getMax();
108  Double_t min = dynamic_cast<RooAbsRealLValue*>(poi)->getMin();
109  Double_t step = (max-min)/(curve_N-1);
110  curve_x=new Double_t[curve_N];
111  for(int i=0; i<curve_N; ++i){
112  curve_x[i]=min+step*i;
113  }
114 // }
115 
116  map<string, std::vector<Double_t> > name_val;
117  for(int i=0; i<curve_N; i++){
118  poi->setVal(curve_x[i]);
119  profile->getVal();
120 
121  TIterator* nuis_params_itr=nuis_params->createIterator();
122  TObject* nuis_params_obj;
123  while((nuis_params_obj=nuis_params_itr->Next())){
124  RooRealVar* nuis_param = dynamic_cast<RooRealVar*>(nuis_params_obj);
125  if(nuis_param) {
126  string name = nuis_param->GetName();
127  if(nuis_params->getSize()==0) continue;
128  if(nuis_param && (! nuis_param->isConstant())){
129  if(name_val.find(name)==name_val.end()) name_val[name]=std::vector<Double_t>(curve_N);
130  name_val[name][i]=nuis_param->getVal();
131 
132  if(i==curve_N-1){
133  TGraph* g = new TGraph(curve_N, curve_x, &(name_val[name].front()));
134  g->SetName((name+"_"+string(poi->GetName())+"_profile").c_str());
135  g->GetXaxis()->SetTitle(poi->GetName());
136  g->GetYaxis()->SetTitle(nuis_param->GetName());
137  g->SetTitle("");
138  list->Add(g);
139  }
140  }
141  }
142  }
143  }
144 
145  delete [] curve_x;
146 
147 
148  delete nll;
149  delete profile;
150  return list;
151 }