Logo ROOT  
Reference Guide
df033_Describe.py File Reference

Namespaces

 df033_Describe
 

Detailed Description

View in nbviewer Open in SWAN Get information about the dataframe with the convenience method Describe.

import ROOT
# Create a dataframe
path = 'root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/Run2012BC_DoubleMuParked_Muons.root'
df = ROOT.RDataFrame('Events', path)
# Describe the state of the dataframe.
# Note that this operation is not running the event loop.
print(df.Describe())
# Build a small analysis studying the invariant mass of dimuon systems.
# See tutorial df102_NanoAODDimuonAnalysis for more information.
df = df.Filter('nMuon == 2')\
.Filter('Muon_charge[0] != Muon_charge[1]')\
.Define('Dimuon_mass', 'InvariantMass(Muon_pt, Muon_eta, Muon_phi, Muon_mass)')\
.Filter('Dimuon_mass > 70')\
.Range(1000)
# Trigger the event loop by asking for the mean of the dimuon mass.
print('\nApproximate mass of the Z boson: {:.2f} GeV\n'.format(
df.Mean('Dimuon_mass').GetValue()))
# Describe again the state of the dataframe.
print(df.Describe())
Dataframe from TChain Events in file root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/Run2012BC_DoubleMuParked_Muons.root
Property Value
-------- -----
Columns in total 6
Columns from defines 0
Event loops run 0
Column Type Origin
------ ---- ------
nMuon UInt_t Dataset
Muon_pt ROOT::VecOps::RVec<Float_t> Dataset
Muon_eta ROOT::VecOps::RVec<Float_t> Dataset
Muon_phi ROOT::VecOps::RVec<Float_t> Dataset
Muon_mass ROOT::VecOps::RVec<Float_t> Dataset
Muon_charge ROOT::VecOps::RVec<Int_t> Dataset
Approximate mass of the Z boson: 91.44 GeV
Dataframe from TChain Events in file root://eospublic.cern.ch//eos/opendata/cms/derived-data/AOD2NanoAODOutreachTool/Run2012BC_DoubleMuParked_Muons.root
Property Value
-------- -----
Columns in total 7
Columns from defines 1
Event loops run 1
Column Type Origin
------ ---- ------
Dimuon_mass float Define
nMuon UInt_t Dataset
Muon_pt ROOT::VecOps::RVec<Float_t> Dataset
Muon_eta ROOT::VecOps::RVec<Float_t> Dataset
Muon_phi ROOT::VecOps::RVec<Float_t> Dataset
Muon_mass ROOT::VecOps::RVec<Float_t> Dataset
Muon_charge ROOT::VecOps::RVec<Int_t> Dataset
Date
March 2021
Author
Stefan Wunsch (KIT, CERN)

Definition in file df033_Describe.py.

Range
Ta Range(0, 0, 1, 1)
ROOT::RDataFrame
ROOT's RDataFrame offers a high level interface for analyses of data stored in TTree,...
Definition: RDataFrame.hxx:42
ROOT::VecOps::Filter
RVec< T > Filter(const RVec< T > &v, F &&f)
Create a new collection with the elements passing the filter expressed by the predicate.
Definition: RVec.hxx:944