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#include <ROOT/RDataFrame.hxx>
17
18#include <TCanvas.h>
19#include <TH1D.h>
20#include <TLatex.h>
21#include <TStyle.h>
22
23#include <cassert>
24#include <cmath>
25#include <iostream>
26#include <memory>
27#include <string>
28#include <vector>
29#include <utility>
30
31constexpr char const *kNTupleFileName = "http://root.cern/files/tutorials/ntpl004_dimuon_v1.root";
32
33using namespace ROOT::VecOps;
34
35void ntpl004_dimuon() {
36 // Use all available CPU cores
38
40
41 // The tutorial is identical to df102_NanoAODDimuonAnalysis except the use of RNTuple.
42
43 // For simplicity, select only events with exactly two muons and require opposite charge
44 auto df_2mu = df.Filter("nMuon == 2", "Events with exactly two muons");
45 auto df_os = df_2mu.Filter("Muon_charge[0] != Muon_charge[1]", "Muons with opposite charge");
46
47 // Compute invariant mass of the dimuon system
48 auto df_mass = df_os.Define("Dimuon_mass", InvariantMass<float>, {"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
49
50 // Make histogram of dimuon mass spectrum
51 auto h = df_mass.Histo1D({"Dimuon_mass", "Dimuon_mass", 30000, 0.25, 300}, "Dimuon_mass");
52
53 // Request cut-flow report
54 auto report = df_mass.Report();
55
56 // Produce plot
58 auto c = new TCanvas("c", "", 800, 700);
59 c->SetLogx(); c->SetLogy();
60
61 h->SetTitle("");
62 h->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)"); h->GetXaxis()->SetTitleSize(0.04);
63 h->GetYaxis()->SetTitle("N_{Events}"); h->GetYaxis()->SetTitleSize(0.04);
64 h->DrawCopy();
65
66 TLatex label; label.SetNDC(true);
67 label.DrawLatex(0.175, 0.740, "#eta");
68 label.DrawLatex(0.205, 0.775, "#rho,#omega");
69 label.DrawLatex(0.270, 0.740, "#phi");
70 label.DrawLatex(0.400, 0.800, "J/#psi");
71 label.DrawLatex(0.415, 0.670, "#psi'");
72 label.DrawLatex(0.485, 0.700, "Y(1,2,3S)");
73 label.DrawLatex(0.755, 0.680, "Z");
74 label.SetTextSize(0.040); label.DrawLatex(0.100, 0.920, "#bf{CMS Open Data}");
75 label.SetTextSize(0.030); label.DrawLatex(0.630, 0.920, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}");
76
77 // Print cut-flow report
78 report->Print();
79}
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
R__EXTERN TStyle * gStyle
Definition TStyle.h:442
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:48
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:49
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:1642
virtual void SetNDC(Bool_t isNDC=kTRUE)
Set NDC mode on if isNDC = kTRUE, off otherwise.
Definition TText.cxx:817
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