Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "RooProfileLL.h"
41
42#include "RooPlot.h"
43
44#include "TCanvas.h"
45#include "TH1F.h"
46#include "TH2F.h"
47#include "TTree.h"
48#include "TMarker.h"
49#include "TStopwatch.h"
50
51#include <iostream>
52
53// use this order for safety on library loading
54using namespace RooFit;
55using namespace RooStats;
56
57void rs401d_FeldmanCousins(bool doFeldmanCousins = false, bool doMCMC = true)
58{
60
61 // to time the macro
62 TStopwatch t;
63 t.Start();
64
65 // Taken from Feldman & Cousins paper, Phys.Rev.D57:3873-3889,1998.
66 // e-Print: physics/9711021 (see page 13.)
67 //
68 // Quantum mechanics dictates that the probability of such a transformation is given by the formula
69 // $P (\nu\mu \rightarrow \nu e ) = sin^2 (2\theta) sin^2 (1.27 \Delta m^2 L /E )$
70 // where P is the probability for a $\nu\mu$ to transform into a $\nu e$ , L is the distance in km between
71 // the creation of the neutrino from meson decay and its interaction in the detector, E is the
72 // neutrino energy in GeV, and $\Delta m^2 = |m^2 - m^2 |$ in $(eV/c^2 )^2$ .
73 //
74 // To demonstrate how this works in practice, and how it compares to alternative approaches
75 // that have been used, we consider a toy model of a typical neutrino oscillation experiment.
76 // The toy model is defined by the following parameters: Mesons are assumed to decay to
77 // neutrinos uniformly in a region 600 m to 1000 m from the detector. The expected background
78 // from conventional $\nu e$ interactions and misidentified $\nu\mu$ interactions is assumed to be 100
79 // events in each of 5 energy bins which span the region from 10 to 60 GeV. We assume that
80 // the $\nu\mu$ flux is such that if $P (\nu\mu \rightarrow \nu e ) = 0.01$ averaged over any bin, then that bin
81 // would
82 // have an expected additional contribution of 100 events due to $\nu\mu \rightarrow \nu e$ oscillations.
83
84 // Make signal model model
85 RooRealVar E("E", "", 15, 10, 60, "GeV");
86 RooRealVar L("L", "", .800, .600, 1.0, "km"); // need these units in formula
87 RooRealVar deltaMSq("deltaMSq", "#Delta m^{2}", 40, 1, 300, "eV/c^{2}");
88 RooRealVar sinSq2theta("sinSq2theta", "sin^{2}(2#theta)", .006, .0, .02);
89 // RooRealVar deltaMSq("deltaMSq","#Delta m^{2}",40,20,70,"eV/c^{2}");
90 // RooRealVar sinSq2theta("sinSq2theta","sin^{2}(2#theta)", .006,.001,.01);
91 // PDF for oscillation only describes deltaMSq dependence, sinSq2theta goes into sigNorm
92 auto oscillationFormula = "std::pow(std::sin(1.27 * x[2] * x[0] / x[1]), 2)";
93 RooGenericPdf PnmuTone("PnmuTone", "P(#nu_{#mu} #rightarrow #nu_{e}", oscillationFormula, {L, E, deltaMSq});
94
95 // only E is observable, so create the signal model by integrating out L
96 RooAbsPdf *sigModel = PnmuTone.createProjection(L);
97
98 // create $ \int dE' dL' P(E',L' | \Delta m^2)$.
99 // Given RooFit will renormalize the PDF in the range of the observables,
100 // the average probability to oscillate in the experiment's acceptance
101 // needs to be incorporated into the extended term in the likelihood.
102 // Do this by creating a RooAbsReal representing the integral and divide by
103 // the area in the E-L plane.
104 // The integral should be over "primed" observables, so we need
105 // an independent copy of PnmuTone not to interfere with the original.
106
107 // Independent copy for Integral
108 RooRealVar EPrime("EPrime", "", 15, 10, 60, "GeV");
109 RooRealVar LPrime("LPrime", "", .800, .600, 1.0, "km"); // need these units in formula
110 RooGenericPdf PnmuTonePrime("PnmuTonePrime", "P(#nu_{#mu} #rightarrow #nu_{e}", oscillationFormula, {LPrime, EPrime, deltaMSq});
112
113 // Getting the flux is a bit tricky. It is more clear to include a cross section term that is not
114 // explicitly referred to in the text, eg.
115 // number events in bin = flux * cross-section for nu_e interaction in E bin * average prob nu_mu osc. to nu_e in bin
116 // let maxEventsInBin = flux * cross-section for nu_e interaction in E bin
117 // maxEventsInBin * 1% chance per bin = 100 events / bin
118 // therefore maxEventsInBin = 10,000.
119 // for 5 bins, this means maxEventsTot = 50,000
120 RooConstVar maxEventsTot("maxEventsTot", "maximum number of sinal events", 50000);
121 RooConstVar inverseArea("inverseArea", "1/(#Delta E #Delta L)",
122 1. / (EPrime.getMax() - EPrime.getMin()) / (LPrime.getMax() - LPrime.getMin()));
123
124 // $sigNorm = maxEventsTot \cdot \int dE dL \frac{P_{oscillate\ in\ experiment}}{Area} \cdot {sin}^2(2\theta)$
126 // bkg = 5 bins * 100 events / bin
127 RooConstVar bkgNorm("bkgNorm", "normalization for background", 500);
128
129 // flat background (0th order polynomial, so no arguments for coefficients)
130 RooPolynomial bkgEShape("bkgEShape", "flat bkg shape", E);
131
132 // total model
134
135 // for debugging, check model tree
136 // model.printCompactTree();
137 // model.graphVizTree("model.dot");
138
139 // turn off some messages
140 RooMsgService::instance().setStreamStatus(0, false);
141 RooMsgService::instance().setStreamStatus(1, false);
142 RooMsgService::instance().setStreamStatus(2, false);
143
144 // --------------------------------------
145 // n events in data to data, simply sum of sig+bkg
146 Int_t nEventsData = bkgNorm.getVal() + sigNorm.getVal();
147 cout << "generate toy data with nEvents = " << nEventsData << endl;
148 // adjust random seed to get a toy dataset similar to one in paper.
149 // Found by trial and error (3 trials, so not very "fine tuned")
150 RooRandom::randomGenerator()->SetSeed(3);
151 // create a toy dataset
152 RooDataSet *data = model.generate(RooArgSet(E), nEventsData);
153
154 // --------------------------------------
155 // make some plots
156 TCanvas *dataCanvas = new TCanvas("dataCanvas");
157 dataCanvas->Divide(2, 2);
158
159 // plot the PDF
160 dataCanvas->cd(1);
161 TH1 *hh = PnmuTone.createHistogram("hh", E, Binning(40), YVar(L, Binning(40)), Scaling(false));
162 hh->SetLineColor(kBlue);
163 hh->SetTitle("True Signal Model");
164 hh->Draw("surf");
165
166 // plot the data with the best fit
167 dataCanvas->cd(2);
168 RooPlot *Eframe = E.frame();
169 data->plotOn(Eframe);
170 model.fitTo(*data, Extended());
171 model.plotOn(Eframe);
172 model.plotOn(Eframe, Components(*sigModel), LineColor(kRed));
173 model.plotOn(Eframe, Components(bkgEShape), LineColor(kGreen));
174 model.plotOn(Eframe);
175 Eframe->SetTitle("toy data with best fit model (and sig+bkg components)");
176 Eframe->Draw();
177
178 // plot the likelihood function
179 dataCanvas->cd(3);
180 std::unique_ptr<RooAbsReal> nll{model.createNLL(*data, Extended(true))};
181 RooProfileLL pll("pll", "", *nll, RooArgSet(deltaMSq, sinSq2theta));
182 // TH1* hhh = nll.createHistogram("hhh",sinSq2theta,Binning(40),YVar(deltaMSq,Binning(40))) ;
183 TH1 *hhh = pll.createHistogram("hhh", sinSq2theta, Binning(40), YVar(deltaMSq, Binning(40)), Scaling(false));
184 hhh->SetLineColor(kBlue);
185 hhh->SetTitle("Likelihood Function");
186 hhh->Draw("surf");
187
188 dataCanvas->Update();
189
190 // --------------------------------------------------------------
191 // show use of Feldman-Cousins utility in RooStats
192 // set the distribution creator, which encodes the test statistic
193 RooArgSet parameters(deltaMSq, sinSq2theta);
194 RooWorkspace *w = new RooWorkspace();
195
197 modelConfig.SetWorkspace(*w);
198 modelConfig.SetPdf(model);
199 modelConfig.SetParametersOfInterest(parameters);
200
202 fc.SetTestSize(.1); // set size of test
203 fc.UseAdaptiveSampling(true);
204 fc.SetNBins(10); // number of points to test per parameter
205
206 // use the Feldman-Cousins tool
209 interval = fc.GetInterval();
210
211 // ---------------------------------------------------------
212 // show use of ProfileLikeihoodCalculator utility in RooStats
214 plc.SetTestSize(.1);
215
216 ConfInterval *plcInterval = plc.GetInterval();
217
218 // --------------------------------------------
219 // show use of MCMCCalculator utility in RooStats
221
222 if (doMCMC) {
223 // turn some messages back on
224 RooMsgService::instance().setStreamStatus(0, true);
225 RooMsgService::instance().setStreamStatus(1, true);
226
228 mcmcWatch.Start();
229
232 mc.SetNumIters(5000);
233 mc.SetNumBurnInSteps(100);
234 mc.SetUseKeys(true);
235 mc.SetTestSize(.1);
236 mc.SetAxes(axisList); // set which is x and y axis in posterior histogram
237 // mc.SetNumBins(50);
238 mcInt = (MCMCInterval *)mc.GetInterval();
239
240 mcmcWatch.Stop();
241 mcmcWatch.Print();
242 }
243 // -------------------------------
244 // make plot of resulting interval
245
246 dataCanvas->cd(4);
247
248 // first plot a small dot for every point tested
249 if (doFeldmanCousins) {
250 RooDataHist *parameterScan = (RooDataHist *)fc.GetPointsToScan();
251 TH2F *hist = parameterScan->createHistogram(deltaMSq,sinSq2theta, 30, 30);
252 // hist->Draw();
253 TH2F *forContour = (TH2F *)hist->Clone();
254
255 // now loop through the points and put a marker if it's in the interval
257 // loop over points to test
258 for (Int_t i = 0; i < parameterScan->numEntries(); ++i) {
259 // get a parameter point from the list of points to test.
260 tmpPoint = (RooArgSet *)parameterScan->get(i)->clone("temp");
261
262 if (interval) {
263 if (interval->IsInInterval(*tmpPoint)) {
264 forContour->SetBinContent(
265 hist->FindBin(tmpPoint->getRealValue("sinSq2theta"), tmpPoint->getRealValue("deltaMSq")), 1);
266 } else {
267 forContour->SetBinContent(
268 hist->FindBin(tmpPoint->getRealValue("sinSq2theta"), tmpPoint->getRealValue("deltaMSq")), 0);
269 }
270 }
271
272 delete tmpPoint;
273 }
274
275 if (interval) {
276 Double_t level = 0.5;
277 forContour->SetContour(1, &level);
278 forContour->SetLineWidth(2);
279 forContour->SetLineColor(kRed);
280 forContour->Draw("cont2,same");
281 }
282 }
283
285 if (mcInt) {
286 cout << "MCMC actual confidence level: " << mcInt->GetActualConfidenceLevel() << endl;
288 mcPlot->SetLineColor(kMagenta);
289 mcPlot->Draw();
290 }
291 dataCanvas->Update();
292
294 plotInt.SetTitle("90% Confidence Intervals");
295 if (mcInt)
296 plotInt.Draw("same");
297 else
298 plotInt.Draw();
299 dataCanvas->Update();
300
301 /// print timing info
302 t.Stop();
303 t.Print();
304}
int Int_t
Signed integer 4 bytes (int)
Definition RtypesCore.h:59
double Double_t
Double 8 bytes.
Definition RtypesCore.h:73
@ kRed
Definition Rtypes.h:67
@ kGreen
Definition Rtypes.h:67
@ kMagenta
Definition Rtypes.h:67
@ kBlue
Definition Rtypes.h:67
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Represents a constant real-valued object.
Definition RooConstVar.h:23
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Container class to hold unbinned data.
Definition RooDataSet.h:32
Implementation of a probability density function that takes a RooArgList of servers and a C++ express...
static RooMsgService & instance()
Return reference to singleton instance.
Plot frame and a container for graphics objects within that frame.
Definition RooPlot.h:43
RooPolynomial implements a polynomial p.d.f of the form.
Represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
Implements the profile likelihood estimator for a given likelihood and set of parameters of interest.
static TRandom * randomGenerator()
Return a pointer to a singleton random-number generator implementation.
Definition RooRandom.cxx:47
Variable that can be changed from the outside.
Definition RooRealVar.h:37
static void enableSilentClipping(bool flag=true)
Enable or disable the silent clipping behavior of RooRealVar::setVal() that was the default in ROOT v...
ConfInterval is an interface class for a generic interval in the RooStats framework.
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.
MCMCInterval is a concrete implementation of the RooStats::ConfInterval interface.
< A class that holds configuration information for a model using a workspace as a store
Definition ModelConfig.h:34
The ProfileLikelihoodCalculator is a concrete implementation of CombinedCalculator (the interface cla...
Persistable container for RooFit projects.
The Canvas class.
Definition TCanvas.h:23
TH1 is the base class of all histogram classes in ROOT.
Definition TH1.h:109
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:3727
TObject * Clone(const char *newname="") const override
Make a complete copy of the underlying object.
Definition TH1.cxx:2786
2-D histogram with a float per channel (see TH1 documentation)
Definition TH2.h:345
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
void Stop()
Stop the stopwatch.
void Print(Option_t *option="") const override
Print the real and cpu time passed between the start and stop events.
RooCmdArg YVar(const RooAbsRealLValue &var, const RooCmdArg &arg={})
RooCmdArg Scaling(bool flag)
RooCmdArg Components(Args_t &&... argsOrArgSet)
RooCmdArg Binning(const RooAbsBinning &binning)
RooCmdArg LineColor(TColorNumber color)
RooArgList L(Args_t &&... args)
Definition RooArgList.h:156
double nll(double pdf, double weight, int binnedL, int doBinOffset)
Definition MathFuncs.h:452
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:72
Namespace for the RooStats classes.
Definition CodegenImpl.h:66
constexpr Double_t E()
Base of natural log: .
Definition TMath.h:96