00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028 #include "TFile.h"
00029 #include "TROOT.h"
00030 #include "TH1F.h"
00031
00032 #include "RooWorkspace.h"
00033 #include "RooAbsData.h"
00034
00035 #include "RooStats/ModelConfig.h"
00036 #include "RooStats/FeldmanCousins.h"
00037 #include "RooStats/ToyMCSampler.h"
00038 #include "RooStats/PointSetInterval.h"
00039 #include "RooStats/ConfidenceBelt.h"
00040
00041
00042 using namespace RooFit;
00043 using namespace RooStats;
00044
00045 void StandardFeldmanCousinsDemo(const char* infile = "",
00046 const char* workspaceName = "proto",
00047 const char* modelConfigName = "ModelConfig",
00048 const char* dataName = "expData"){
00049
00050
00051
00052
00053
00054 const char* filename = "";
00055 if (!strcmp(infile,""))
00056 filename = "results/example_channel1_GammaExample_model.root";
00057 else
00058 filename = infile;
00059
00060 TFile *file = TFile::Open(filename);
00061
00062
00063 if(!file && strcmp(infile,"")){
00064 cout <<"file not found" << endl;
00065 return;
00066 }
00067
00068
00069 if(!file ){
00070
00071 cout <<"will run standard hist2workspace example"<<endl;
00072 gROOT->ProcessLine(".! prepareHistFactory .");
00073 gROOT->ProcessLine(".! hist2workspace config/example.xml");
00074 cout <<"\n\n---------------------"<<endl;
00075 cout <<"Done creating example input"<<endl;
00076 cout <<"---------------------\n\n"<<endl;
00077 }
00078
00079
00080 file = TFile::Open(filename);
00081 if(!file){
00082
00083 cout << "Not able to run hist2workspace to create example input" <<endl;
00084 return;
00085 }
00086
00087
00088
00089
00090
00091
00092
00093 RooWorkspace* w = (RooWorkspace*) file->Get(workspaceName);
00094 if(!w){
00095 cout <<"workspace not found" << endl;
00096 return;
00097 }
00098
00099
00100 ModelConfig* mc = (ModelConfig*) w->obj(modelConfigName);
00101
00102
00103 RooAbsData* data = w->data(dataName);
00104
00105
00106 if(!data || !mc){
00107 w->Print();
00108 cout << "data or ModelConfig was not found" <<endl;
00109 return;
00110 }
00111
00112
00113
00114
00115
00116
00117 FeldmanCousins fc(*data,*mc);
00118 fc.SetConfidenceLevel(0.95);
00119
00120 fc.UseAdaptiveSampling(true);
00121 fc.SetNBins(10);
00122 fc.CreateConfBelt(true);
00123
00124
00125
00126
00127
00128
00129
00130
00131 if(!mc->GetPdf()->canBeExtended()){
00132 if(data->numEntries()==1)
00133 fc.FluctuateNumDataEntries(false);
00134 else
00135 cout <<"Not sure what to do about this model" <<endl;
00136 }
00137
00138
00139
00140
00141
00142
00143
00144
00145 PointSetInterval* interval = fc.GetInterval();
00146 ConfidenceBelt* belt = fc.GetConfidenceBelt();
00147
00148
00149 RooRealVar* firstPOI = (RooRealVar*) mc->GetParametersOfInterest()->first();
00150 cout << "\n95% interval on " <<firstPOI->GetName()<<" is : ["<<
00151 interval->LowerLimit(*firstPOI) << ", "<<
00152 interval->UpperLimit(*firstPOI) <<"] "<<endl;
00153
00154
00155
00156
00157
00158 RooDataSet* parameterScan = (RooDataSet*) fc.GetPointsToScan();
00159 RooArgSet* tmpPoint;
00160
00161
00162 TH1F* histOfThresholds = new TH1F("histOfThresholds","",
00163 parameterScan->numEntries(),
00164 firstPOI->getMin(),
00165 firstPOI->getMax());
00166
00167
00168
00169
00170 for(Int_t i=0; i<parameterScan->numEntries(); ++i){
00171 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
00172 double arMax = belt->GetAcceptanceRegionMax(*tmpPoint);
00173 double arMin = belt->GetAcceptanceRegionMax(*tmpPoint);
00174 double poiVal = tmpPoint->getRealValue(firstPOI->GetName()) ;
00175 histOfThresholds->Fill(poiVal,arMax);
00176 }
00177 histOfThresholds->SetMinimum(0);
00178 histOfThresholds->Draw();
00179
00180 }