Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
df037_TTreeEventMatching.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_dataframe
3/// \notebook -nodraw
4///
5/// This example shows processing of a TTree-based dataset with horizontal
6/// concatenations (friends) and event matching (based on TTreeIndex). In case
7/// the current event being processed does not match one (or more) of the friend
8/// datasets, one can use the FilterAvailable and DefaultValueFor functionalities
9/// to act upon the situation.
10///
11/// \macro_code
12/// \macro_output
13///
14/// \date September 2024
15/// \author Vincenzo Eduardo Padulano (CERN)
16#include <TChain.h>
17#include <TFile.h>
18#include <TTree.h>
19#include <TTreeIndex.h>
20
21#include <ROOT/RDataFrame.hxx>
22
23#include <iostream>
24#include <numeric>
25
26// A helper class to create the dataset for the tutorial below.
27struct Dataset {
28
29 constexpr static auto fMainFile{"df037_TTreeEventMatching_C_main.root"};
30 constexpr static auto fAuxFile1{"df037_TTreeEventMatching_C_aux_1.root"};
31 constexpr static auto fAuxFile2{"df037_TTreeEventMatching_C_aux_2.root"};
32 constexpr static auto fMainTreeName{"events"};
33 constexpr static auto fAuxTreeName1{"auxdata_1"};
34 constexpr static auto fAuxTreeName2{"auxdata_2"};
35
36 Dataset()
37 {
38 {
39 TFile f(fMainFile, "RECREATE");
40 TTree mainTree(fMainTreeName, fMainTreeName);
41 int idx;
42 int x;
43 mainTree.Branch("idx", &idx, "idx/I");
44 mainTree.Branch("x", &x, "x/I");
45
46 idx = 1;
47 x = 1;
48 mainTree.Fill();
49 idx = 2;
50 x = 2;
51 mainTree.Fill();
52 idx = 3;
53 x = 3;
54 mainTree.Fill();
55
56 mainTree.Write();
57 }
58 {
59 // The first auxiliary file has matching indices 1 and 2, but not 3
60 TFile f(fAuxFile1, "RECREATE");
61 TTree auxTree(fAuxTreeName1, fAuxTreeName1);
62 int y;
63 int idx;
64 auxTree.Branch("idx", &idx, "idx/I");
65 auxTree.Branch("y", &y, "y/I");
66
67 idx = 1;
68 y = 4;
69 auxTree.Fill();
70 idx = 2;
71 y = 5;
72 auxTree.Fill();
73
74 auxTree.Write();
75 }
76 {
77 // The second auxiliary file has matching indices 1 and 3, but not 2
78 TFile f(fAuxFile2, "RECREATE");
79 TTree auxTree(fAuxTreeName2, fAuxTreeName2);
80 int z;
81 int idx;
82 auxTree.Branch("idx", &idx, "idx/I");
83 auxTree.Branch("z", &z, "z/I");
84
85 idx = 1;
86 z = 6;
87 auxTree.Fill();
88 idx = 3;
89 z = 7;
90 auxTree.Fill();
91
92 auxTree.Write();
93 }
94 }
95
96 ~Dataset()
97 {
98 std::remove(fMainFile);
99 std::remove(fAuxFile1);
100 std::remove(fAuxFile2);
101 }
102};
103
105{
106 // Create the dataset: one main TTree and two auxiliary. The 'idx' branch
107 // is used as the index to match events between the trees.
108 // - The main tree has 3 entries, with 'idx' values (1, 2, 3).
109 // - The first auxiliary tree has 2 entries, with 'idx' values (1, 2).
110 // - The second auxiliary tree has 2 entries, with 'idx' values (1, 3).
111 // The two auxiliary trees are concatenated horizontally with the main one.
112 Dataset dataset{};
113 TChain mainChain{dataset.fMainTreeName};
114 mainChain.Add(dataset.fMainFile);
115
116 TChain auxChain1(dataset.fAuxTreeName1);
117 auxChain1.Add(dataset.fAuxFile1);
118 auxChain1.BuildIndex("idx");
119
120 TChain auxChain2(dataset.fAuxTreeName2);
121 auxChain2.Add(dataset.fAuxFile2);
122 auxChain2.BuildIndex("idx");
123
124 mainChain.AddFriend(&auxChain1);
125 mainChain.AddFriend(&auxChain2);
126
127 // Create an RDataFrame to process the input dataset. The DefaultValueFor and
128 // FilterAvailable functionalities can be used to decide what to do for
129 // the events that do not match entirely according to the index column 'idx'
130 ROOT::RDataFrame df{mainChain};
131
132 const std::string auxTree1ColIdx = std::string(dataset.fAuxTreeName1) + ".idx";
133 const std::string auxTree1ColY = std::string(dataset.fAuxTreeName1) + ".y";
134 const std::string auxTree2ColIdx = std::string(dataset.fAuxTreeName2) + ".idx";
135 const std::string auxTree2ColZ = std::string(dataset.fAuxTreeName2) + ".z";
136
137 constexpr static auto defaultValue = std::numeric_limits<int>::min();
138
139 // Example 1: provide default values for all columns in case there was no
140 // match
141 auto display1 = df.DefaultValueFor(auxTree1ColIdx, defaultValue)
142 .DefaultValueFor(auxTree1ColY, defaultValue)
143 .DefaultValueFor(auxTree2ColIdx, defaultValue)
144 .DefaultValueFor(auxTree2ColZ, defaultValue)
145 .Display<int, int, int, int, int, int>(
146 {"idx", auxTree1ColIdx, auxTree2ColIdx, "x", auxTree1ColY, auxTree2ColZ});
147
148 // Example 2: skip the entire entry when there was no match for a column
149 // in the first auxiliary tree, but keep the entries when there is no match
150 // in the second auxiliary tree and provide a default value for those
151 auto display2 = df.DefaultValueFor(auxTree2ColIdx, defaultValue)
152 .DefaultValueFor(auxTree2ColZ, defaultValue)
153 .FilterAvailable(auxTree1ColY)
154 .Display<int, int, int, int, int, int>(
155 {"idx", auxTree1ColIdx, auxTree2ColIdx, "x", auxTree1ColY, auxTree2ColZ});
156
157 // Example 3: Keep entries from the main tree for which there is no
158 // corresponding match in entries of the first auxiliary tree
159 auto display3 = df.FilterMissing(auxTree1ColIdx).Display<int, int>({"idx", "x"});
160
161 std::cout << "Example 1: provide default values for all columns\n";
162 display1->Print();
163 std::cout << "Example 2: skip the entry only when the first auxiliary tree does not match\n";
164 display2->Print();
165 std::cout << "Example 3: keep entries from the main tree for which there is no match in the auxiliary tree\n";
166 display3->Print();
167}
#define f(i)
Definition RSha256.hxx:104
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
A chain is a collection of files containing TTree objects.
Definition TChain.h:33
A ROOT file is an on-disk file, usually with extension .root, that stores objects in a file-system-li...
Definition TFile.h:53
A TTree represents a columnar dataset.
Definition TTree.h:79
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17