Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rs401d_FeldmanCousins.py
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 (C++ version), and P. P. (Python translation)
16
17import ROOT
18
19
20def rs401d_FeldmanCousins(doFeldmanCousins=False, doMCMC=True):
21
23
24 # to time the macro
25 t = ROOT.TStopwatch()
26 t.Start()
27
28 # Taken from Feldman & Cousins paper, Phys.Rev.D57:3873-3889,1998.
29 # e-Print: physics/9711021 (see page 13.)
30 #
31 # Quantum mechanics dictates that the probability of such a transformation is given by the formula
32 # $P (\nu\mu \rightarrow \nu e ) = sin^2 (2\theta) sin^2 (1.27 \Delta m^2 L /E )$
33 # where P is the probability for a $\nu\mu$ to transform into a $\nu e$ , L is the distance in km between
34 # the creation of the neutrino from meson decay and its interaction in the detector, E is the
35 # neutrino energy in GeV, and $\Delta m^2 = |m^2 - m^2 |$ in $(eV/c^2 )^2$ .
36 #
37 # To demonstrate how this works in practice, and how it compares to alternative approaches
38 # that have been used, we consider a toy model of a typical neutrino oscillation experiment.
39 # The toy model is defined by the following parameters: Mesons are assumed to decay to
40 # neutrinos uniformly in a region 600 m to 1000 m from the detector. The expected background
41 # from conventional $\nu e$ interactions and misidentified $\nu\mu$ interactions is assumed to be 100
42 # events in each of 5 energy bins which span the region from 10 to 60 GeV. We assume that
43 # the $\nu\mu$ flux is such that if $P (\nu\mu \rightarrow \nu e ) = 0.01$ averaged over any bin, then that bin
44 # would
45 # have an expected additional contribution of 100 events due to $\nu\mu \rightarrow \nu e$ oscillations.
46
47 # Make signal model model
48 E = ROOT.RooRealVar("E", "", 15, 10, 60, "GeV")
49 L = ROOT.RooRealVar("L", "", 0.800, 0.600, 1.0, "km") # need these units in formula
50 deltaMSq = ROOT.RooRealVar("deltaMSq", "#Delta m^{2}", 40, 1, 300, "eV/c^{2}")
51 sinSq2theta = ROOT.RooRealVar("sinSq2theta", "sin^{2}(2#theta)", 0.006, 0.0, 0.02)
52 # RooRealVar deltaMSq("deltaMSq","#Delta m^{2}",40,20,70,"eV/c^{2}");
53 # RooRealVar sinSq2theta("sinSq2theta","sin^{2}(2#theta)", .006,.001,.01);
54 # PDF for oscillation only describes deltaMSq dependence, sinSq2theta goes into sigNorm
55 oscillationFormula = "std::pow(std::sin(1.27 * x[2] * x[0] / x[1]), 2)"
56 PnmuTone = ROOT.RooGenericPdf("PnmuTone", "P(#nu_{#mu} #rightarrow #nu_{e}", oscillationFormula, [L, E, deltaMSq])
57
58 # only E is observable, so create the signal model by integrating out L
59 sigModel = PnmuTone.createProjection(L)
60
61 # create $ \int dE' dL' P(E',L' | \Delta m^2)$.
62 # Given RooFit will renormalize the PDF in the range of the observables,
63 # the average probability to oscillate in the experiment's acceptance
64 # needs to be incorporated into the extended term in the likelihood.
65 # Do this by creating a RooAbsReal representing the integral and divide by
66 # the area in the E-L plane.
67 # The integral should be over "primed" observables, so we need
68 # an independent copy of PnmuTone not to interfere with the original.
69
70 # Independent copy for Integral
71 EPrime = ROOT.RooRealVar("EPrime", "", 15, 10, 60, "GeV")
72 LPrime = ROOT.RooRealVar("LPrime", "", 0.800, 0.600, 1.0, "km") # need these units in formula
73 PnmuTonePrime = ROOT.RooGenericPdf(
74 "PnmuTonePrime", "P(#nu_{#mu} #rightarrow #nu_{e}", oscillationFormula, [LPrime, EPrime, deltaMSq]
75 )
76 intProbToOscInExp = PnmuTonePrime.createIntegral([EPrime, LPrime])
77
78 # Getting the flux is a bit tricky. It is more clear to include a cross section term that is not
79 # explicitly referred to in the text, eg.
80 # number events in bin = flux * cross-section for nu_e interaction in E bin * average prob nu_mu osc. to nu_e in bin
81 # let maxEventsInBin = flux * cross-section for nu_e interaction in E bin
82 # maxEventsInBin * 1% chance per bin = 100 events / bin
83 # therefore maxEventsInBin = 10,000.
84 # for 5 bins, this means maxEventsTot = 50,000
85 maxEventsTot = ROOT.RooConstVar("maxEventsTot", "maximum number of sinal events", 50000)
86 inverseArea = ROOT.RooConstVar(
87 "inverseArea",
88 "1/(#Delta E #Delta L)",
90 )
91
92 # $sigNorm = maxEventsTot \cdot \int dE dL \frac{P_{oscillate\ in\ experiment}}{Area} \cdot {sin}^2(2\theta)$
93 sigNorm = ROOT.RooProduct("sigNorm", "", [maxEventsTot, intProbToOscInExp, inverseArea, sinSq2theta])
94 # bkg = 5 bins * 100 events / bin
95 bkgNorm = ROOT.RooConstVar("bkgNorm", "normalization for background", 500)
96
97 # flat background (0th order polynomial, so no arguments for coefficients)
98 bkgEShape = ROOT.RooPolynomial("bkgEShape", "flat bkg shape", E)
99
100 # total model
101 model = ROOT.RooAddPdf("model", "", [sigModel, bkgEShape], [sigNorm, bkgNorm])
102
103 # for debugging, check model tree
104 # model.printCompactTree();
105 # model.graphVizTree("model.dot");
106
107 # turn off some messages
108 ROOT.RooMsgService.instance().setStreamStatus(0, False)
109 ROOT.RooMsgService.instance().setStreamStatus(1, False)
110 ROOT.RooMsgService.instance().setStreamStatus(2, False)
111
112 # --------------------------------------
113 # n events in data to data, simply sum of sig+bkg
114 nEventsData = bkgNorm.getVal() + sigNorm.getVal()
115 print("generate toy data with nEvents = ", nEventsData)
116 # adjust random seed to get a toy dataset similar to one in paper.
117 # Found by trial and error (3 trials, so not very "fine tuned")
119 # create a toy dataset
120 data = model.generate(E, nEventsData)
121
122 # --------------------------------------
123 # make some plots
124 dataCanvas = ROOT.TCanvas("dataCanvas")
126
127 # plot the PDF
129 hh = PnmuTone.createHistogram("hh", E, ROOT.RooFit.YVar(L, ROOT.RooFit.Binning(40)), Binning=40, Scaling=False)
131 hh.SetTitle("True Signal Model")
132 hh.Draw("surf")
135 # dataCanvas.SaveAs("rs.1.png")
136 # plot the data with the best fit
138 Eframe = E.frame()
139 data.plotOn(Eframe)
140 model.fitTo(data, Extended=True)
141 model.plotOn(Eframe)
142 model.plotOn(Eframe, Components=sigModel, LineColor="r")
143 model.plotOn(Eframe, Components=bkgEShape, LineColor="g")
144 model.plotOn(Eframe)
145 Eframe.SetTitle("toy data with best fit model (and sig+bkg components)")
149 # dataCanvas.SaveAs("rs.2.png")
150
151 # plot the likelihood function
153 nll = model.createNLL(data, Extended=True)
154 pll = ROOT.RooProfileLL("pll", "", nll, [deltaMSq, sinSq2theta])
155 # hhh = nll.createHistogram("hhh",sinSq2theta, Binning(40), YVar(deltaMSq,Binning(40)))
157 "hhh", sinSq2theta, ROOT.RooFit.YVar(deltaMSq, ROOT.RooFit.Binning(40)), Binning=40, Scaling=False
158 )
160 hhh.SetTitle("Likelihood Function")
161 hhh.Draw("surf")
162
165 dataCanvas.SaveAs("3.png")
166
167 # --------------------------------------------------------------
168 # show use of Feldman-Cousins utility in RooStats
169 # set the distribution creator, which encodes the test statistic
170 parameters = ROOT.RooArgSet(deltaMSq, sinSq2theta)
172
173 modelConfig = ROOT.RooFit.ModelConfig()
175 modelConfig.SetPdf(model)
177
178 fc = ROOT.RooStats.FeldmanCousins(data, modelConfig)
179 fc.SetTestSize(0.1) # set size of test
181 fc.SetNBins(10) # number of points to test per parameter
182
183 # use the Feldman-Cousins tool
184 interval = 0
185 if doFeldmanCousins:
186 interval = fc.GetInterval()
187
188 # ---------------------------------------------------------
189 # show use of ProfileLikeihoodCalculator utility in RooStats
190 plc = ROOT.RooStats.ProfileLikelihoodCalculator(data, modelConfig)
191 plc.SetTestSize(0.1)
192
193 plcInterval = plc.GetInterval()
194
195 # --------------------------------------------
196 # show use of MCMCCalculator utility in RooStats
197 mcInt = ROOT.kNone
198
199 if doMCMC:
200 # turn some messages back on
201 ROOT.RooMsgService.instance().setStreamStatus(0, True)
202 ROOT.RooMsgService.instance().setStreamStatus(1, True)
203
204 mcmcWatch = ROOT.TStopwatch()
206
207 axisList = ROOT.RooArgList(deltaMSq, sinSq2theta)
208 mc = ROOT.RooStats.MCMCCalculator(data, modelConfig)
209 mc.SetNumIters(5000)
211 mc.SetUseKeys(True)
212 mc.SetTestSize(0.1)
213 mc.SetAxes(axisList) # set which is x and y axis in posterior histogram
214 # mc.SetNumBins(50);
215 mcInt = mc.GetInterval()
216
219
220 # -------------------------------
221 # make plot of resulting interval
222
224
225 # first plot a small dot for every point tested
226 if doFeldmanCousins:
227 parameterScan = fc.GetPointsToScan()
228 hist = parameterScan.createHistogram(deltaMSq, sinSq2theta, 30, 30)
229 # hist.Draw()
230 forContour = hist.Clone()
231
232 # now loop through the points and put a marker if it's in the interval
233 tmpPoint = RooArgSet()
234 # loop over points to test
236 # get a parameter point from the list of points to test.
237 tmpPoint = parameterScan.get(i).clone("temp")
238
239 if interval:
240 if interval.IsInInterval(tmpPoint):
242 hist.FindBin(tmpPoint.getRealValue("sinSq2theta"), tmpPoint.getRealValue("deltaMSq")), 1
243 )
244 else:
246 hist.FindBin(tmpPoint.getRealValue("sinSq2theta"), tmpPoint.getRealValue("deltaMSq")), 0
247 )
248
249 del tmpPoint
250
251 if interval:
252 level = 0.5
253 forContour.SetContour(1, level)
256 forContour.Draw("cont2,same")
257
258 mcPlot = ROOT.kNone
259 if mcInt:
260 print(f"MCMC actual confidence level: ", mcInt.GetActualConfidenceLevel())
261 mcPlot = ROOT.RooStats.MCMCIntervalPlot(mcInt)
264
266
267 plotInt = ROOT.RooStats.LikelihoodIntervalPlot(plcInterval)
268
269 plotInt.SetTitle("90% Confidence Intervals")
270 if mcInt:
271 plotInt.Draw("same")
272 else:
274
276
277 dataCanvas.SaveAs("rs401d_FeldmanCousins.1.pdf")
278
279 # print timing info
280 t.Stop()
281 t.Print()
282
283
284rs401d_FeldmanCousins(doFeldmanCousins=False, doMCMC=True)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24