Logo ROOT   6.16/01
Reference Guide
rs401d_FeldmanCousins.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_roostats
3/// \notebook
4/// 'Neutrino Oscillation Example from Feldman & Cousins'
5///
6/// This tutorial shows a more complex example using the FeldmanCousins utility
7/// to create a confidence interval for a toy neutrino oscillation experiment.
8/// The example attempts to faithfully reproduce the toy example described in Feldman & Cousins'
9/// original paper, Phys.Rev.D57:3873-3889,1998.
10///
11/// \macro_image
12/// \macro_output
13/// \macro_code
14///
15/// \author Kyle Cranmer
16
17#include "RooGlobalFunc.h"
26
27#include "RooDataSet.h"
28#include "RooDataHist.h"
29#include "RooRealVar.h"
30#include "RooConstVar.h"
31#include "RooAddition.h"
32#include "RooProduct.h"
33#include "RooProdPdf.h"
34#include "RooAddPdf.h"
35
36#include "TROOT.h"
37#include "RooPolynomial.h"
38#include "RooRandom.h"
39
40#include "RooNLLVar.h"
41#include "RooProfileLL.h"
42
43#include "RooPlot.h"
44
45#include "TCanvas.h"
46#include "TH1F.h"
47#include "TH2F.h"
48#include "TTree.h"
49#include "TMarker.h"
50#include "TStopwatch.h"
51
52#include <iostream>
53
54// PDF class created for this macro
55#if !defined(__CINT__) || defined(__MAKECINT__)
56#include "../tutorials/roostats/NuMuToNuE_Oscillation.h"
57#include "../tutorials/roostats/NuMuToNuE_Oscillation.cxx" // so that it can be executed directly
58#else
59#include "../tutorials/roostats/NuMuToNuE_Oscillation.cxx+" // so that it can be executed directly
60#endif
61
62// use this order for safety on library loading
63using namespace RooFit;
64using namespace RooStats ;
65
66
67void rs401d_FeldmanCousins(bool doFeldmanCousins=false, bool doMCMC = true)
68{
69
70 // to time the macro
71 TStopwatch t;
72 t.Start();
73
74
75
76 // Taken from Feldman & Cousins paper, Phys.Rev.D57:3873-3889,1998.
77 // e-Print: physics/9711021 (see page 13.)
78 //
79 // Quantum mechanics dictates that the probability of such a transformation is given by the formula
80 // $P (\nu\mu \rightarrow \nu e ) = sin^2 (2\theta) sin^2 (1.27 \Delta m^2 L /E )$
81 // where P is the probability for a $\nu\mu$ to transform into a $\nu e$ , L is the distance in km between
82 // the creation of the neutrino from meson decay and its interaction in the detector, E is the
83 // neutrino energy in GeV, and $\Delta m^2 = |m^2 - m^2 |$ in $(eV/c^2 )^2$ .
84 //
85 // To demonstrate how this works in practice, and how it compares to alternative approaches
86 // that have been used, we consider a toy model of a typical neutrino oscillation experiment.
87 // The toy model is defined by the following parameters: Mesons are assumed to decay to
88 // neutrinos uniformly in a region 600 m to 1000 m from the detector. The expected background
89 // from conventional $\nu e$ interactions and misidentified $\nu\mu$ interactions is assumed to be 100
90 // events in each of 5 energy bins which span the region from 10 to 60 GeV. We assume that
91 // the $\nu\mu$ flux is such that if $P (\nu\mu \rightarrow \nu e ) = 0.01$ averaged over any bin, then that bin would
92 // have an expected additional contribution of 100 events due to $\nu\mu \rightarrow \nu e$ oscillations.
93
94
95 // Make signal model model
96 RooRealVar E("E","", 15,10,60,"GeV");
97 RooRealVar L("L","", .800,.600, 1.0,"km"); // need these units in formula
98 RooRealVar deltaMSq("deltaMSq","#Delta m^{2}",40,1,300,"eV/c^{2}");
99 RooRealVar sinSq2theta("sinSq2theta","sin^{2}(2#theta)", .006,.0,.02);
100 //RooRealVar deltaMSq("deltaMSq","#Delta m^{2}",40,20,70,"eV/c^{2}");
101 // RooRealVar sinSq2theta("sinSq2theta","sin^{2}(2#theta)", .006,.001,.01);
102 // PDF for oscillation only describes deltaMSq dependence, sinSq2theta goes into sigNorm
103 // 1) The code for this PDF was created by issuing these commands
104 // root [0] RooClassFactory x
105 // root [1] x.makePdf("NuMuToNuE_Oscillation","L,E,deltaMSq","","pow(sin(1.27*deltaMSq*L/E),2)")
106 NuMuToNuE_Oscillation PnmuTone("PnmuTone","P(#nu_{#mu} #rightarrow #nu_{e}",L,E,deltaMSq);
107
108 // only E is observable, so create the signal model by integrating out L
109 RooAbsPdf* sigModel = PnmuTone.createProjection(L);
110
111 // create $ \int dE' dL' P(E',L' | \Delta m^2)$.
112 // Given RooFit will renormalize the PDF in the range of the observables,
113 // the average probability to oscillate in the experiment's acceptance
114 // needs to be incorporated into the extended term in the likelihood.
115 // Do this by creating a RooAbsReal representing the integral and divide by
116 // the area in the E-L plane.
117 // The integral should be over "primed" observables, so we need
118 // an independent copy of PnmuTone not to interfere with the original.
119
120 // Independent copy for Integral
121 RooRealVar EPrime("EPrime","", 15,10,60,"GeV");
122 RooRealVar LPrime("LPrime","", .800,.600, 1.0,"km"); // need these units in formula
123 NuMuToNuE_Oscillation PnmuTonePrime("PnmuTonePrime","P(#nu_{#mu} #rightarrow #nu_{e}",
124 LPrime,EPrime,deltaMSq);
125 RooAbsReal* intProbToOscInExp = PnmuTonePrime.createIntegral(RooArgSet(EPrime,LPrime));
126
127 // Getting the flux is a bit tricky. It is more clear to include a cross section term that is not
128 // explicitly referred to in the text, eg.
129 // number events in bin = flux * cross-section for nu_e interaction in E bin * average prob nu_mu osc. to nu_e in bin
130 // let maxEventsInBin = flux * cross-section for nu_e interaction in E bin
131 // maxEventsInBin * 1% chance per bin = 100 events / bin
132 // therefore maxEventsInBin = 10,000.
133 // for 5 bins, this means maxEventsTot = 50,000
134 RooConstVar maxEventsTot("maxEventsTot","maximum number of sinal events",50000);
135 RooConstVar inverseArea("inverseArea","1/(#Delta E #Delta L)",
136 1./(EPrime.getMax()-EPrime.getMin())/(LPrime.getMax()-LPrime.getMin()));
137
138 // $sigNorm = maxEventsTot \cdot \int dE dL \frac{P_{oscillate\ in\ experiment}}{Area} \cdot {sin}^2(2\theta)$
139 RooProduct sigNorm("sigNorm", "", RooArgSet(maxEventsTot, *intProbToOscInExp, inverseArea, sinSq2theta));
140 // bkg = 5 bins * 100 events / bin
141 RooConstVar bkgNorm("bkgNorm","normalization for background",500);
142
143 // flat background (0th order polynomial, so no arguments for coefficients)
144 RooPolynomial bkgEShape("bkgEShape","flat bkg shape", E);
145
146 // total model
147 RooAddPdf model("model","",RooArgList(*sigModel,bkgEShape),
148 RooArgList(sigNorm,bkgNorm));
149
150 // for debugging, check model tree
151 // model.printCompactTree();
152 // model.graphVizTree("model.dot");
153
154
155 // turn off some messages
159
160
161 // --------------------------------------
162 // n events in data to data, simply sum of sig+bkg
163 Int_t nEventsData = bkgNorm.getVal()+sigNorm.getVal();
164 cout << "generate toy data with nEvents = " << nEventsData << endl;
165 // adjust random seed to get a toy dataset similar to one in paper.
166 // Found by trial and error (3 trials, so not very "fine tuned")
168 // create a toy dataset
169 RooDataSet* data = model.generate(RooArgSet(E), nEventsData);
170
171 // --------------------------------------
172 // make some plots
173 TCanvas* dataCanvas = new TCanvas("dataCanvas");
174 dataCanvas->Divide(2,2);
175
176 // plot the PDF
177 dataCanvas->cd(1);
178 TH1* hh = PnmuTone.createHistogram("hh",E,Binning(40),YVar(L,Binning(40)),Scaling(kFALSE)) ;
179 hh->SetLineColor(kBlue) ;
180 hh->SetTitle("True Signal Model");
181 hh->Draw("surf");
182
183 // plot the data with the best fit
184 dataCanvas->cd(2);
185 RooPlot* Eframe = E.frame();
186 data->plotOn(Eframe);
187 model.fitTo(*data, Extended());
188 model.plotOn(Eframe);
189 model.plotOn(Eframe,Components(*sigModel),LineColor(kRed));
190 model.plotOn(Eframe,Components(bkgEShape),LineColor(kGreen));
191 model.plotOn(Eframe);
192 Eframe->SetTitle("toy data with best fit model (and sig+bkg components)");
193 Eframe->Draw();
194
195 // plot the likelihood function
196 dataCanvas->cd(3);
197 RooNLLVar nll("nll", "nll", model, *data, Extended());
198 RooProfileLL pll("pll", "", nll, RooArgSet(deltaMSq, sinSq2theta));
199 // TH1* hhh = nll.createHistogram("hhh",sinSq2theta,Binning(40),YVar(deltaMSq,Binning(40))) ;
200 TH1* hhh = pll.createHistogram("hhh",sinSq2theta,Binning(40),YVar(deltaMSq,Binning(40)),Scaling(kFALSE)) ;
201 hhh->SetLineColor(kBlue) ;
202 hhh->SetTitle("Likelihood Function");
203 hhh->Draw("surf");
204
205 dataCanvas->Update();
206
207
208
209 // --------------------------------------------------------------
210 // show use of Feldman-Cousins utility in RooStats
211 // set the distribution creator, which encodes the test statistic
212 RooArgSet parameters(deltaMSq, sinSq2theta);
213 RooWorkspace* w = new RooWorkspace();
214
215 ModelConfig modelConfig;
216 modelConfig.SetWorkspace(*w);
217 modelConfig.SetPdf(model);
218 modelConfig.SetParametersOfInterest(parameters);
219
220 RooStats::FeldmanCousins fc(*data, modelConfig);
221 fc.SetTestSize(.1); // set size of test
222 fc.UseAdaptiveSampling(true);
223 fc.SetNBins(10); // number of points to test per parameter
224
225 // use the Feldman-Cousins tool
226 ConfInterval* interval = 0;
227 if(doFeldmanCousins)
228 interval = fc.GetInterval();
229
230
231 // ---------------------------------------------------------
232 // show use of ProfileLikeihoodCalculator utility in RooStats
234 plc.SetTestSize(.1);
235
236 ConfInterval* plcInterval = plc.GetInterval();
237
238 // --------------------------------------------
239 // show use of MCMCCalculator utility in RooStats
240 MCMCInterval* mcInt = NULL;
241
242 if (doMCMC) {
243 // turn some messages back on
246
247 TStopwatch mcmcWatch;
248 mcmcWatch.Start();
249
250 RooArgList axisList(deltaMSq, sinSq2theta);
251 MCMCCalculator mc(*data, modelConfig);
252 mc.SetNumIters(5000);
253 mc.SetNumBurnInSteps(100);
254 mc.SetUseKeys(true);
255 mc.SetTestSize(.1);
256 mc.SetAxes(axisList); // set which is x and y axis in posterior histogram
257 //mc.SetNumBins(50);
258 mcInt = (MCMCInterval*)mc.GetInterval();
259
260 mcmcWatch.Stop();
261 mcmcWatch.Print();
262 }
263 // -------------------------------
264 // make plot of resulting interval
265
266 dataCanvas->cd(4);
267
268 // first plot a small dot for every point tested
269 if (doFeldmanCousins) {
270 RooDataHist* parameterScan = (RooDataHist*) fc.GetPointsToScan();
271 TH2F* hist = (TH2F*) parameterScan->createHistogram("sinSq2theta:deltaMSq",30,30);
272 // hist->Draw();
273 TH2F* forContour = (TH2F*)hist->Clone();
274
275 // now loop through the points and put a marker if it's in the interval
276 RooArgSet* tmpPoint;
277 // loop over points to test
278 for(Int_t i=0; i<parameterScan->numEntries(); ++i){
279 // get a parameter point from the list of points to test.
280 tmpPoint = (RooArgSet*) parameterScan->get(i)->clone("temp");
281
282 if (interval){
283 if (interval->IsInInterval( *tmpPoint ) ) {
284 forContour->SetBinContent( hist->FindBin(tmpPoint->getRealValue("sinSq2theta"),
285 tmpPoint->getRealValue("deltaMSq")), 1);
286 }else{
287 forContour->SetBinContent( hist->FindBin(tmpPoint->getRealValue("sinSq2theta"),
288 tmpPoint->getRealValue("deltaMSq")), 0);
289 }
290 }
291
292
293 delete tmpPoint;
294 }
295
296 if (interval){
297 Double_t level=0.5;
298 forContour->SetContour(1,&level);
299 forContour->SetLineWidth(2);
300 forContour->SetLineColor(kRed);
301 forContour->Draw("cont2,same");
302 }
303 }
304
305 MCMCIntervalPlot* mcPlot = NULL;
306 if (mcInt) {
307 cout << "MCMC actual confidence level: "
308 << mcInt->GetActualConfidenceLevel() << endl;
309 mcPlot = new MCMCIntervalPlot(*mcInt);
310 mcPlot->SetLineColor(kMagenta);
311 mcPlot->Draw();
312 }
313 dataCanvas->Update();
314
315 LikelihoodIntervalPlot plotInt((LikelihoodInterval*)plcInterval);
316 plotInt.SetTitle("90% Confidence Intervals");
317 if (mcInt)
318 plotInt.Draw("same");
319 else
320 plotInt.Draw();
321 dataCanvas->Update();
322
323 /// print timing info
324 t.Stop();
325 t.Print();
326
327
328}
329
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: RtypesCore.h:88
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
@ kRed
Definition: Rtypes.h:63
@ kGreen
Definition: Rtypes.h:63
@ kMagenta
Definition: Rtypes.h:63
@ kBlue
Definition: Rtypes.h:63
static struct mg_connection * fc(struct mg_context *ctx)
Definition: civetweb.c:3728
TH1 * createHistogram(const char *name, const RooAbsRealLValue &xvar, const RooCmdArg &arg1=RooCmdArg::none(), const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Calls createHistogram(const char *name, const RooAbsRealLValue& xvar, const RooLinkedList& argList) c...
Definition: RooAbsData.cxx:613
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
virtual RooAbsPdf * createProjection(const RooArgSet &iset)
Return a p.d.f that represent a projection of this p.d.f integrated over given observables.
Definition: RooAbsPdf.cxx:2914
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
RooAbsReal * createIntegral(const RooArgSet &iset, const RooCmdArg &arg1, const RooCmdArg &arg2=RooCmdArg::none(), const RooCmdArg &arg3=RooCmdArg::none(), const RooCmdArg &arg4=RooCmdArg::none(), const RooCmdArg &arg5=RooCmdArg::none(), const RooCmdArg &arg6=RooCmdArg::none(), const RooCmdArg &arg7=RooCmdArg::none(), const RooCmdArg &arg8=RooCmdArg::none()) const
Create an object that represents the integral of the function over one or more observables listed in ...
Definition: RooAbsReal.cxx:502
RooAddPdf is an efficient implementation of a sum of PDFs of the form.
Definition: RooAddPdf.h:29
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
Double_t getRealValue(const char *name, Double_t defVal=0, Bool_t verbose=kFALSE) const
Get value of a RooAbsReal stored in set with given name.
Definition: RooArgSet.cxx:472
virtual TObject * clone(const char *newname) const
Definition: RooArgSet.h:84
RooConstVar represent a constant real-valued object.
Definition: RooConstVar.h:25
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
virtual Int_t numEntries() const
Return the number of bins.
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
RooDataSet is a container class to hold unbinned data.
Definition: RooDataSet.h:31
void setStreamStatus(Int_t id, Bool_t active)
(De)Activate stream with given unique ID
static RooMsgService & instance()
Return reference to singleton instance.
Class RooNLLVar implements a a -log(likelihood) calculation from a dataset and a PDF.
Definition: RooNLLVar.h:26
A RooPlot is a plot frame and a container for graphics objects within that frame.
Definition: RooPlot.h:41
void SetTitle(const char *name)
Set the title of the RooPlot to 'title'.
Definition: RooPlot.cxx:1104
virtual void Draw(Option_t *options=0)
Draw this plot and all of the elements it contains.
Definition: RooPlot.cxx:558
RooPolynomial implements a polynomial p.d.f of the form.
Definition: RooPolynomial.h:28
RooProduct a RooAbsReal implementation that represent the product of a given set of other RooAbsReal ...
Definition: RooProduct.h:32
Class RooProfileLL implements the profile likelihood estimator for a given likelihood and set of para...
Definition: RooProfileLL.h:26
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition: RooRandom.cxx:54
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
ConfInterval is an interface class for a generic interval in the RooStats framework.
Definition: ConfInterval.h:35
virtual Bool_t IsInInterval(const RooArgSet &) const =0
check if given point is in the interval
The FeldmanCousins class (like the Feldman-Cousins technique) is essentially a specific configuration...
This class provides simple and straightforward utilities to plot a LikelihoodInterval object.
LikelihoodInterval is a concrete implementation of the RooStats::ConfInterval interface.
Bayesian Calculator estimating an interval or a credible region using the Markov-Chain Monte Carlo me...
This class provides simple and straightforward utilities to plot a MCMCInterval object.
void SetLineColor(Color_t color)
void Draw(const Option_t *options=NULL)
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
Definition: MCMCInterval.h:30
virtual Double_t GetActualConfidenceLevel()
virtual Double_t GetKeysPdfCutoff() { return fKeysCutoff; }
ModelConfig is a simple class that holds configuration information specifying how a model should be u...
Definition: ModelConfig.h:30
virtual void SetWorkspace(RooWorkspace &ws)
Definition: ModelConfig.h:66
virtual void SetParametersOfInterest(const RooArgSet &set)
Definition: ModelConfig.h:99
virtual void SetPdf(const RooAbsPdf &pdf)
Set the Pdf, add to the the workspace if not already there.
Definition: ModelConfig.h:81
ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface class f...
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition: TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition: TAttLine.h:40
The Canvas class.
Definition: TCanvas.h:31
virtual void Update()
Update canvas pad buffers.
Definition: TCanvas.cxx:2286
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:693
The TH1 histogram class.
Definition: TH1.h:56
virtual void SetTitle(const char *title)
See GetStatOverflows for more information.
Definition: TH1.cxx:6217
TObject * Clone(const char *newname=0) const
Make a complete copy of the underlying object.
Definition: TH1.cxx:2657
virtual void SetContour(Int_t nlevels, const Double_t *levels=0)
Set the number and values of contour levels.
Definition: TH1.cxx:7813
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2974
virtual Int_t FindBin(Double_t x, Double_t y=0, Double_t z=0)
Return Global bin number corresponding to x,y,z.
Definition: TH1.cxx:3572
2-D histogram with a float per channel (see TH1 documentation)}
Definition: TH2.h:250
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content.
Definition: TH2.cxx:2500
virtual void Divide(Int_t nx=1, Int_t ny=1, Float_t xmargin=0.01, Float_t ymargin=0.01, Int_t color=0)
Automatic pad generation by division.
Definition: TPad.cxx:1162
virtual void SetSeed(ULong_t seed=0)
Set the random generator seed.
Definition: TRandom.cxx:589
Stopwatch class.
Definition: TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void Stop()
Stop the stopwatch.
Definition: TStopwatch.cxx:77
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
RooCmdArg Scaling(Bool_t flag)
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg=RooCmdArg::none())
RooCmdArg Extended(Bool_t flag=kTRUE)
RooCmdArg Components(const RooArgSet &compSet)
RooCmdArg LineColor(Color_t color)
@(#)root/roostats:$Id$ Author: George Lewis, Kyle Cranmer
Definition: Asimov.h:20
static constexpr double L
constexpr Double_t E()
Base of natural log:
Definition: TMath.h:97