Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
ntpl004_dimuon.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Mini-Analysis on CMS OpenData with RDataFrame.

This tutorial illustrates that analyzing data with RDataFrame works the same for both TTree data and RNTuple data. The RNTuple data are converted from the Events tree in http://root.cern/files/NanoAOD_DoubleMuon_CMS2011OpenData.root Based on RDataFrame's df102_NanoAODDimuonAnalysis.C

#include <TCanvas.h>
#include <TH1D.h>
#include <TLatex.h>
#include <TStyle.h>
#include <cassert>
#include <cmath>
#include <iostream>
#include <memory>
#include <string>
#include <vector>
#include <utility>
constexpr char const *kNTupleFileName = "http://root.cern/files/tutorials/ntpl004_dimuon_v1.root";
using namespace ROOT::VecOps;
void ntpl004_dimuon() {
// Use all available CPU cores
ROOT::RDataFrame df("Events", kNTupleFileName);
// The tutorial is identical to df102_NanoAODDimuonAnalysis except the use of RNTuple.
// For simplicity, select only events with exactly two muons and require opposite charge
auto df_2mu = df.Filter("nMuon == 2", "Events with exactly two muons");
auto df_os = df_2mu.Filter("Muon_charge[0] != Muon_charge[1]", "Muons with opposite charge");
// Compute invariant mass of the dimuon system
auto df_mass = df_os.Define("Dimuon_mass", InvariantMass<float>, {"Muon_pt", "Muon_eta", "Muon_phi", "Muon_mass"});
// Make histogram of dimuon mass spectrum
auto h = df_mass.Histo1D({"Dimuon_mass", "Dimuon_mass", 30000, 0.25, 300}, "Dimuon_mass");
// Request cut-flow report
auto report = df_mass.Report();
// Produce plot
gStyle->SetOptStat(0); gStyle->SetTextFont(42);
auto c = new TCanvas("c", "", 800, 700);
c->SetLogx(); c->SetLogy();
h->SetTitle("");
h->GetXaxis()->SetTitle("m_{#mu#mu} (GeV)"); h->GetXaxis()->SetTitleSize(0.04);
h->GetYaxis()->SetTitle("N_{Events}"); h->GetYaxis()->SetTitleSize(0.04);
h->DrawCopy();
TLatex label; label.SetNDC(true);
label.DrawLatex(0.175, 0.740, "#eta");
label.DrawLatex(0.205, 0.775, "#rho,#omega");
label.DrawLatex(0.270, 0.740, "#phi");
label.DrawLatex(0.400, 0.800, "J/#psi");
label.DrawLatex(0.415, 0.670, "#psi'");
label.DrawLatex(0.485, 0.700, "Y(1,2,3S)");
label.DrawLatex(0.755, 0.680, "Z");
label.SetTextSize(0.040); label.DrawLatex(0.100, 0.920, "#bf{CMS Open Data}");
label.SetTextSize(0.030); label.DrawLatex(0.630, 0.920, "#sqrt{s} = 8 TeV, L_{int} = 11.6 fb^{-1}");
// Print cut-flow report
report->Print();
}
#define c(i)
Definition RSha256.hxx:101
#define h(i)
Definition RSha256.hxx:106
externTStyle * gStyle
Definition TStyle.h:442
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
virtual void SetTextSize(Float_t tsize=1)
Set the text size.
Definition TAttText.h:53
The Canvas class.
Definition TCanvas.h:23
TLatex * DrawLatex(Double_t x, Double_t y, const char *text)
virtual void SetNDC(Bool_t isNDC=kTRUE)
void Print(Option_t *option="") const override
Print TNamed name and title.
Common_t InvariantMass(const RVec< T0 > &pt, const RVec< T1 > &eta, const RVec< T2 > &phi, const RVec< T3 > &mass)
Return the invariant mass of multiple particles given the collections of the quantities transverse mo...
Definition RVec.hxx:3164
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:613
Date
April 2019
Author
The ROOT Team

Definition in file ntpl004_dimuon.C.