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

Detailed Description

View in nbviewer Open in SWAN
Show how NanoAOD files can be processed with RDataFrame.

This tutorial illustrates how NanoAOD files can be processed with ROOT dataframes. The NanoAOD-like input files are filled with 66 mio. events from CMS OpenData containing muon candidates part of 2012 dataset (DOI: 10.7483/OPENDATA.CMS.YLIC.86ZZ and DOI: 10.7483/OPENDATA.CMS.M5AD.Y3V3). The macro matches muon pairs and produces an histogram of the dimuon mass spectrum showing resonances up to the Z mass. Note that the bump at 30 GeV is not a resonance but a trigger effect.

More details about the dataset can be found on the CERN Open Data portal.

#include "ROOT/RVec.hxx"
#include "TCanvas.h"
#include "TH1D.h"
#include "TLatex.h"
#include "TStyle.h"
using namespace ROOT::VecOps;
void df102_NanoAODDimuonAnalysis()
{
// Enable multi-threading
// Create dataframe from NanoAOD files
ROOT::RDataFrame df("Events", "root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/"
"Run2012BC_DoubleMuParked_Muons.root");
// Add ProgressBar
// 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. Note how we can set title and axis labels in one go
auto h = df_mass.Histo1D({"Dimuon_mass", "Dimuon mass;m_{#mu#mu} (GeV);N_{Events}", 30000, 0.25, 300}, "Dimuon_mass");
// Request cut-flow report
auto report = df.Report();
// Produce plot
gStyle->SetOptStat(0); gStyle->SetTextFont(42);
auto c = new TCanvas("c", "", 800, 700);
c->SetLogx(); c->SetLogy();
h->GetXaxis()->SetTitleSize(0.04);
h->GetYaxis()->SetTitleSize(0.04);
h->DrawClone();
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}");
c->SaveAs("dimuon_spectrum.pdf");
// Print cut-flow report
report->Print();
}
int main()
{
df102_NanoAODDimuonAnalysis();
}
int main()
Definition Prototype.cxx:12
#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)
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 AddProgressBar(ROOT::RDF::RNode df)
Add ProgressBar to a ROOT::RDF::RNode.
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
|=> | [Elapsed: 0:10m files: 1 / 1 events: 13195000 / 61540413 1.32e+06 evt/s remaining ca.: 0:36m]
|=======> | [Elapsed: 0:20m files: 1 / 1 events: 43979000 / 61540413 2.20e+06 evt/s remaining ca.: 0:07m]
[Total elapsed time: 0:26m processed files: 1 processed events: 61540413 2.37e+06 evt/s]
Events with exactly two muons: pass=31104343 all=61540413 -- eff=50.54 % cumulative eff=50.54 %
Muons with opposite charge: pass=24067843 all=31104343 -- eff=77.38 % cumulative eff=39.11 %
Date
August 2018
Author
Stefan Wunsch (KIT, CERN)

Definition in file df102_NanoAODDimuonAnalysis.C.