00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016 #include "RooGlobalFunc.h"
00017 #include "RooStats/ConfInterval.h"
00018 #include "RooStats/FeldmanCousins.h"
00019 #include "RooStats/ProfileLikelihoodCalculator.h"
00020 #include "RooStats/MCMCCalculator.h"
00021 #include "RooStats/UniformProposal.h"
00022 #include "RooStats/LikelihoodIntervalPlot.h"
00023 #include "RooStats/MCMCIntervalPlot.h"
00024 #include "RooStats/MCMCInterval.h"
00025
00026 #include "RooDataSet.h"
00027 #include "RooDataHist.h"
00028 #include "RooRealVar.h"
00029 #include "RooConstVar.h"
00030 #include "RooAddition.h"
00031 #include "RooProduct.h"
00032 #include "RooProdPdf.h"
00033 #include "RooAddPdf.h"
00034
00035 #include "TROOT.h"
00036 #include "RooPolynomial.h"
00037 #include "RooRandom.h"
00038
00039 #include "RooNLLVar.h"
00040 #include "RooProfileLL.h"
00041
00042 #include "RooPlot.h"
00043
00044 #include "TCanvas.h"
00045 #include "TH1F.h"
00046 #include "TH2F.h"
00047 #include "TTree.h"
00048 #include "TMarker.h"
00049 #include "TStopwatch.h"
00050
00051 #include <iostream>
00052
00053
00054 #if !defined(__CINT__) || defined(__MAKECINT__)
00055 #include "../tutorials/roostats/NuMuToNuE_Oscillation.h"
00056 #include "../tutorials/roostats/NuMuToNuE_Oscillation.cxx"
00057 #else
00058 #include "../tutorials/roostats/NuMuToNuE_Oscillation.cxx+"
00059 #endif
00060
00061
00062 using namespace RooFit ;
00063 using namespace RooStats ;
00064
00065
00066 void rs401d_FeldmanCousins(bool doFeldmanCousins=false, bool doMCMC = true)
00067 {
00068
00069
00070 TStopwatch t;
00071 t.Start();
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083
00084
00085
00086
00087
00088
00089
00090
00091
00092
00093
00094
00095 RooRealVar E("E","", 15,10,60,"GeV");
00096 RooRealVar L("L","", .800,.600, 1.0,"km");
00097 RooRealVar deltaMSq("deltaMSq","#Delta m^{2}",40,1,300,"eV/c^{2}");
00098 RooRealVar sinSq2theta("sinSq2theta","sin^{2}(2#theta)", .006,.0,.02);
00099
00100
00101
00102
00103
00104
00105 NuMuToNuE_Oscillation PnmuTone("PnmuTone","P(#nu_{#mu} #rightarrow #nu_{e}",L,E,deltaMSq);
00106
00107
00108 RooAbsPdf* sigModel = PnmuTone.createProjection(L);
00109
00110
00111
00112
00113
00114
00115
00116
00117
00118
00119
00120 RooRealVar EPrime("EPrime","", 15,10,60,"GeV");
00121 RooRealVar LPrime("LPrime","", .800,.600, 1.0,"km");
00122 NuMuToNuE_Oscillation PnmuTonePrime("PnmuTonePrime","P(#nu_{#mu} #rightarrow #nu_{e}",
00123 LPrime,EPrime,deltaMSq);
00124 RooAbsReal* intProbToOscInExp = PnmuTonePrime.createIntegral(RooArgSet(EPrime,LPrime));
00125
00126
00127
00128
00129
00130
00131
00132
00133 RooConstVar maxEventsTot("maxEventsTot","maximum number of sinal events",50000);
00134 RooConstVar inverseArea("inverseArea","1/(#Delta E #Delta L)",
00135 1./(EPrime.getMax()-EPrime.getMin())/(LPrime.getMax()-LPrime.getMin()));
00136
00137
00138 RooProduct sigNorm("sigNorm", "", RooArgSet(maxEventsTot, *intProbToOscInExp, inverseArea, sinSq2theta));
00139
00140 RooConstVar bkgNorm("bkgNorm","normalization for background",500);
00141
00142
00143 RooPolynomial bkgEShape("bkgEShape","flat bkg shape", E);
00144
00145
00146 RooAddPdf model("model","",RooArgList(*sigModel,bkgEShape),
00147 RooArgList(sigNorm,bkgNorm));
00148
00149
00150
00151
00152
00153
00154
00155 RooMsgService::instance().setStreamStatus(0,kFALSE);
00156 RooMsgService::instance().setStreamStatus(1,kFALSE);
00157 RooMsgService::instance().setStreamStatus(2,kFALSE);
00158
00159
00160
00161
00162 Int_t nEventsData = bkgNorm.getVal()+sigNorm.getVal();
00163 cout << "generate toy data with nEvents = " << nEventsData << endl;
00164
00165
00166 RooRandom::randomGenerator()->SetSeed(3);
00167
00168 RooDataSet* data = model.generate(RooArgSet(E), nEventsData);
00169
00170
00171
00172 TCanvas* dataCanvas = new TCanvas("dataCanvas");
00173 dataCanvas->Divide(2,2);
00174
00175
00176 dataCanvas->cd(1);
00177 TH1* hh = PnmuTone.createHistogram("hh",E,Binning(40),YVar(L,Binning(40)),Scaling(kFALSE)) ;
00178 hh->SetLineColor(kBlue) ;
00179 hh->SetTitle("True Signal Model");
00180 hh->Draw("surf");
00181
00182
00183 dataCanvas->cd(2);
00184 RooPlot* Eframe = E.frame();
00185 data->plotOn(Eframe);
00186 model.fitTo(*data, Extended());
00187 model.plotOn(Eframe);
00188 model.plotOn(Eframe,Components(*sigModel),LineColor(kRed));
00189 model.plotOn(Eframe,Components(bkgEShape),LineColor(kGreen));
00190 model.plotOn(Eframe);
00191 Eframe->SetTitle("toy data with best fit model (and sig+bkg components)");
00192 Eframe->Draw();
00193
00194
00195 dataCanvas->cd(3);
00196 RooNLLVar nll("nll", "nll", model, *data, Extended());
00197 RooProfileLL pll("pll", "", nll, RooArgSet(deltaMSq, sinSq2theta));
00198
00199 TH1* hhh = pll.createHistogram("hhh",sinSq2theta,Binning(40),YVar(deltaMSq,Binning(40)),Scaling(kFALSE)) ;
00200 hhh->SetLineColor(kBlue) ;
00201 hhh->SetTitle("Likelihood Function");
00202 hhh->Draw("surf");
00203
00204 dataCanvas->Update();
00205
00206
00207
00208
00209
00210
00211 RooArgSet parameters(deltaMSq, sinSq2theta);
00212 RooWorkspace* w = new RooWorkspace();
00213
00214 ModelConfig modelConfig;
00215 modelConfig.SetWorkspace(*w);
00216 modelConfig.SetPdf(model);
00217 modelConfig.SetParametersOfInterest(parameters);
00218
00219 RooStats::FeldmanCousins fc(*data, modelConfig);
00220 fc.SetTestSize(.1);
00221 fc.UseAdaptiveSampling(true);
00222 fc.SetNBins(10);
00223
00224
00225 ConfInterval* interval = 0;
00226 if(doFeldmanCousins)
00227 interval = fc.GetInterval();
00228
00229
00230
00231
00232 RooStats::ProfileLikelihoodCalculator plc(*data, modelConfig);
00233 plc.SetTestSize(.1);
00234
00235 ConfInterval* plcInterval = plc.GetInterval();
00236
00237
00238
00239 MCMCInterval* mcInt = NULL;
00240
00241 if (doMCMC) {
00242
00243 RooMsgService::instance().setStreamStatus(0,kTRUE);
00244 RooMsgService::instance().setStreamStatus(1,kTRUE);
00245
00246 TStopwatch mcmcWatch;
00247 mcmcWatch.Start();
00248
00249 RooArgList axisList(deltaMSq, sinSq2theta);
00250 MCMCCalculator mc(*data, modelConfig);
00251 mc.SetNumIters(5000);
00252 mc.SetNumBurnInSteps(100);
00253 mc.SetUseKeys(true);
00254 mc.SetTestSize(.1);
00255 mc.SetAxes(axisList);
00256
00257 mcInt = (MCMCInterval*)mc.GetInterval();
00258
00259 mcmcWatch.Stop();
00260 mcmcWatch.Print();
00261 }
00262
00263
00264
00265 dataCanvas->cd(4);
00266
00267
00268 if (doFeldmanCousins) {
00269 RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
00270 TH2F* hist = (TH2F*) parameterScan->createHistogram("sinSq2theta:deltaMSq",30,30);
00271
00272 TH2F* forContour = (TH2F*)hist->Clone();
00273
00274
00275 RooArgSet* tmpPoint;
00276
00277 for(Int_t i=0; i<parameterScan->numEntries(); ++i){
00278
00279 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
00280
00281 if (interval){
00282 if (interval->IsInInterval( *tmpPoint ) ) {
00283 forContour->SetBinContent( hist->FindBin(tmpPoint->getRealValue("sinSq2theta"),
00284 tmpPoint->getRealValue("deltaMSq")), 1);
00285 }else{
00286 forContour->SetBinContent( hist->FindBin(tmpPoint->getRealValue("sinSq2theta"),
00287 tmpPoint->getRealValue("deltaMSq")), 0);
00288 }
00289 }
00290
00291
00292 delete tmpPoint;
00293 }
00294
00295 if (interval){
00296 Double_t level=0.5;
00297 forContour->SetContour(1,&level);
00298 forContour->SetLineWidth(2);
00299 forContour->SetLineColor(kRed);
00300 forContour->Draw("cont2,same");
00301 }
00302 }
00303
00304 MCMCIntervalPlot* mcPlot = NULL;
00305 if (mcInt) {
00306 cout << "MCMC actual confidence level: "
00307 << mcInt->GetActualConfidenceLevel() << endl;
00308 mcPlot = new MCMCIntervalPlot(*mcInt);
00309 mcPlot->SetLineColor(kMagenta);
00310 mcPlot->Draw();
00311 }
00312 dataCanvas->Update();
00313
00314 LikelihoodIntervalPlot plotInt((LikelihoodInterval*)plcInterval);
00315 plotInt.SetTitle("90% Confidence Intervals");
00316 if (mcInt)
00317 plotInt.Draw("same");
00318 else
00319 plotInt.Draw();
00320 dataCanvas->Update();
00321
00322
00323 t.Stop();
00324 t.Print();
00325
00326
00327 }
00328