Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
ntpl004_dimuon.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ntuple
3/// \notebook
4/// Mini-Analysis on CMS OpenData with RDataFrame.
5/// This tutorial illustrates that analyzing data with RDataFrame works the same
6/// for both TTree data and RNTuple data. The RNTuple data are converted from the Events tree
7/// in http://root.cern/files/NanoAOD_DoubleMuon_CMS2011OpenData.root
8/// Based on RDataFrame's df102_NanoAODDimuonAnalysis.C
9///
10/// \macro_image
11/// \macro_code
12///
13/// \date April 2019
14/// \author The ROOT Team
15
16// NOTE: The RNTuple classes are experimental at this point.
17// Functionality, interface, and data format is still subject to changes.
18// Do not use for real data!
19
20#include <ROOT/RDataFrame.hxx>
21#include <ROOT/RNTupleDS.hxx>
22#include <ROOT/RVec.hxx>
23
24#include <TCanvas.h>
25#include <TH1D.h>
26#include <TLatex.h>
27#include <TStyle.h>
28
29#include <cassert>
30#include <cmath>
31#include <iostream>
32#include <memory>
33#include <string>
34#include <vector>
35#include <utility>
36
37// Import classes from experimental namespace for the time being
38using RNTupleDS = ROOT::Experimental::RNTupleDS;
39
40constexpr char const* kNTupleFileName = "http://root.cern/files/tutorials/ntpl004_dimuon_v1rc2.root";
41
42using namespace ROOT::VecOps;
43
44void ntpl004_dimuon() {
45 // Use all available CPU cores
47
48 ROOT::RDataFrame df("Events", kNTupleFileName);
49
50 // The tutorial is identical to df102_NanoAODDimuonAnalysis except the use of RNTuple.
51
52 // For simplicity, select only events with exactly two muons and require opposite charge
53 auto df_2mu = df.Filter("nMuon == 2", "Events with exactly two muons");
54 auto df_os = df_2mu.Filter("Muon_charge[0] != Muon_charge[1]", "Muons with opposite charge");
55
56 // Compute invariant mass of the dimuon system
57 auto df_mass = df_os.Define("Dimuon_mass", InvariantMass<float>, {"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
58
59 // Make histogram of dimuon mass spectrum
60 auto h = df_mass.Histo1D({"Dimuon_mass", "Dimuon_mass", 30000, 0.25, 300}, "Dimuon_mass");
61
62 // Request cut-flow report
63 auto report = df_mass.Report();
64
65 // Produce plot
67 auto c = new TCanvas("c", "", 800, 700);
68 c->SetLogx(); c->SetLogy();
69
70 h->SetTitle("");
71 h->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)"); h->GetXaxis()->SetTitleSize(0.04);
72 h->GetYaxis()->SetTitle("N_{Events}"); h->GetYaxis()->SetTitleSize(0.04);
73 h->DrawCopy();
74
75 TLatex label; label.SetNDC(true);
76 label.DrawLatex(0.175, 0.740, "#eta");
77 label.DrawLatex(0.205, 0.775, "#rho,#omega");
78 label.DrawLatex(0.270, 0.740, "#phi");
79 label.DrawLatex(0.400, 0.800, "J/#psi");
80 label.DrawLatex(0.415, 0.670, "#psi'");
81 label.DrawLatex(0.485, 0.700, "Y(1,2,3S)");
82 label.DrawLatex(0.755, 0.680, "Z");
83 label.SetTextSize(0.040); label.DrawLatex(0.100, 0.920, "#bf{CMS Open Data}");
84 label.SetTextSize(0.030); label.DrawLatex(0.630, 0.920, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}");
85
86 // Print cut-flow report
87 report->Print();
88}
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
The RDataSource implementation for RNTuple.
Definition RNTupleDS.hxx:43
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
virtual void SetTextFont(Font_t tfont=62)
Set the text font.
Definition TAttText.h:46
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:47
The Canvas class.
Definition TCanvas.h:23
To draw Mathematical Formula.
Definition TLatex.h:18
TLatex * DrawLatex(Double_t x, Double_t y, const char *text)
Make a copy of this object with the new parameters And copy object attributes.
Definition TLatex.cxx:1943
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1636
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition TText.cxx:823
void Print(Option_t *option="") const override
Dump this text with its attributes.
Definition TText.cxx:788
void EnableImplicitMT(UInt_t numthreads=0)
Enable ROOT's implicit multi-threading for all objects and methods that provide an internal paralleli...
Definition TROOT.cxx:539