Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMVA_SOFIE_GNN_Application.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_ml
3/// \notebook -nodraw
4/// Macro evaluating a GNN model which was generated with the Parser macro
5///
6/// \macro_code
7///
8/// \author
9
10// need to add include path to find generated model file
11#ifdef __CLING__
13#endif
14
15#include "encoder.hxx"
16#include "core.hxx"
17#include "decoder.hxx"
18#include "output_transform.hxx"
19
20#include "TMVA/SOFIE_common.hxx"
21#include "TRandom3.h"
22#include "TH1.h"
23#include "TCanvas.h"
24#include "TFile.h"
25#include "TTree.h"
26#include "TSystem.h"
27#include "ROOT/RDataFrame.hxx"
28
29using namespace TMVA::Experimental;
30using namespace TMVA::Experimental::SOFIE;
31
32double check_mem(std::string s = ""){
34 printf("%s - ",s.c_str());
36 printf(" Rmem = %8.3f MB, Vmem = %8.f3 MB \n",
37 p.fMemResident /1024., /// convert memory from kB to MB
38 p.fMemVirtual /1024.
39 );
40 return p.fMemResident / 1024.;
41}
42
43
44template<class T>
45void PrintTensor(RTensor<T> & t) {
46 std::cout << " shape : " << ConvertShapeToString(t.GetShape()) << " size : " << t.GetSize() << "\n";
47 auto & shape = t.GetShape();
48 auto p = t.GetData();
49 size_t nrows = (shape.size() > 1) ? shape[0] : 1;
50 size_t ncols = (shape.size() > 1) ? t.GetStrides()[0] : shape[0];
51 for (size_t i = 0; i < nrows; i++) {
52 for (size_t j = 0; j < ncols; j++) {
53 if (j==ncols-1) {
54 if (j>10) std::cout << "... ";
55 std::cout << *p << std::endl;
56 }
57 else if (j<10)
58 std::cout << *p << ", ";
59 p++;
60 }
61 }
62 std::cout << std::endl;
63}
64void Print(GNN_Data & d, std::string txt = "") {
65 if (!txt.empty()) std::cout << std::endl << txt << std::endl;
66 std::cout << "node data:"; PrintTensor(d.node_data);
67 std::cout << "edge data:"; PrintTensor(d.edge_data);
68 std::cout << "global data:"; PrintTensor(d.global_data);
69 std::cout << "edge index:"; PrintTensor(d.edge_index);
70}
71
72struct SOFIE_GNN {
73 bool verbose = false;
74 TMVA_SOFIE_encoder::Session encoder;
75 TMVA_SOFIE_core::Session core;
76 TMVA_SOFIE_decoder::Session decoder;
77 TMVA_SOFIE_output_transform::Session output_transform;
78
79 std::vector<GNN_Data> Infer(const GNN_Data & data, int nsteps) {
80 // infer function
81 auto input_data = Copy(data);
82 if (verbose) Print(input_data,"input_data");
83 encoder.infer(input_data);
84 // latent0 is result of encoder. Need to copy because this stays the same
85 auto latent0 = Copy(input_data);
86 GNN_Data latent = input_data; // this can be a view
87 std::vector<GNN_Data> outputData;
88 for (int i = 0; i < nsteps; i++) {
89 if (verbose) Print(latent0,"input decoded data");
90 if (verbose) Print(latent,"latent data");
92 if (verbose) Print(core_input, "after concatenate");
93 core.infer(core_input);
94 if (verbose) Print(core_input, "after core inference");
95 // here I need to copy
97 decoder.infer(core_input);
99 outputData.push_back(Copy(core_input));
100 }
101 return outputData;
102 }
103
104 SOFIE_GNN(bool v = false) : verbose(v) {}
105
106};
107
108const int num_max_nodes = 10;
109const int num_max_edges = 30;
110const int NODE_FEATURE_SIZE = 4;
111const int EDGE_FEATURE_SIZE = 4;
112const int GLOBAL_FEATURE_SIZE = 1;
113
114std::vector<GNN_Data> GenerateData(int nevts, int seed) {
115 TRandom3 r(seed);
116 std::vector<GNN_Data> dataSet;
117 dataSet.reserve(nevts);
118 for (int i = 0; i < nevts; i++) {
119 // generate first number of nodes and edges
120 // size_t num_nodes = num_max_nodes;//r.Integer(num_max_nodes-2) + 2;
121 // size_t num_edges = num_max_edges; //r.Integer(num_max_edges-1) + 1;
122 size_t num_nodes = r.Integer(num_max_nodes-2) + 2;
123 size_t num_edges = r.Integer(num_max_edges-1) + 1;
124 GNN_Data gd;
125 gd.node_data = RTensor<float>({num_nodes, NODE_FEATURE_SIZE});
126 gd.edge_data = RTensor<float>({num_edges, EDGE_FEATURE_SIZE});
127 gd.global_data = RTensor<float>({1, GLOBAL_FEATURE_SIZE});
128 gd.edge_index = RTensor<int>({2, num_edges});
129
130 auto genValue = [&]() { return r.Rndm()*10 -5; };
131 auto genLink = [&] () { return r.Integer(num_nodes);};
132 std::generate(gd.node_data.begin(), gd.node_data.end(), genValue);
133 std::generate(gd.edge_data.begin(), gd.edge_data.end(), genValue);
134 std::generate(gd.global_data.begin(), gd.global_data.end(), genValue);
135 std::generate(gd.edge_index.begin(), gd.edge_index.end(), genLink);
136 dataSet.emplace_back(gd);
137 }
138 return dataSet;
139}
140
141std::vector<GNN_Data> ReadData(std::string treename, std::string filename) {
142 bool verbose = false;
149 int nevts = ndata.GetPtr()->size();
150 std::vector<GNN_Data> dataSet;
151 dataSet.reserve(nevts);
152 for (int i = 0; i < nevts; i++) {
153 GNN_Data gd;
154 auto & n = (*(ndata.GetPtr()))[i];
155 size_t num_nodes = n.size()/NODE_FEATURE_SIZE;
156 auto & e = (*(edata.GetPtr()))[i];
157 size_t num_edges = e.size()/EDGE_FEATURE_SIZE;
158 auto & g = (*(gdata.GetPtr()))[i];
159 gd.node_data = RTensor<float>(n.data(), {num_nodes, NODE_FEATURE_SIZE});
160 gd.edge_data = RTensor<float>(e.data(), {num_edges, EDGE_FEATURE_SIZE});
161 gd.global_data = RTensor<float>(g.data(), {1, GLOBAL_FEATURE_SIZE});
162 gd.edge_index = RTensor<int>({2, num_edges});
163 auto & r = (*(rdata.GetPtr()))[i];
164 auto & s = (*(sdata.GetPtr()))[i];
165 // need to copy receivers/senders in edge_index tensor
166 std::copy(r.begin(), r.end(), gd.edge_index.GetData());
167 std::copy(s.begin(), s.end(), gd.edge_index.GetData()+num_edges);
168
169 dataSet.emplace_back(Copy(gd)); // need to copy data in vector to own
170 if (i < 1 && verbose) Print(dataSet[i],"Input for Event" + std::to_string(i));
171 }
172 return dataSet;
173}
174
175
176void TMVA_SOFIE_GNN_Application (bool verbose = false)
177{
178 check_mem("Initial memory");
180 check_mem("After creating GNN");
181
182
183 const int seed = 111;
184 const int nproc_steps = 5;
185 // generate the input data
186
187 int nevts;
188 //std::cout << "generating data\n";
189 //nevts = 100;
190 //auto inputData = GenerateData(nevts, seed);
191
192 std::cout << "reading data\n";
193 auto inputData = ReadData("gdata","graph_data.root");
194 nevts = inputData.size();
195
196 //std::cout << "padding data\n";
197 //PadData(inputData) ;
198
199 auto h1 = new TH1D("h1","SOFIE Node data",40,1,0);
200 auto h2 = new TH1D("h2","SOFIE Edge data",40,1,0);
201 auto h3 = new TH1D("h3","SOFIE Global data",40,1,0);
202 std::cout << "doing inference...\n";
203
204
205 check_mem("Before evaluating");
206 TStopwatch w; w.Start();
207 for (int i = 0; i < nevts; i++) {
208 auto result = gnn.Infer(inputData[i], nproc_steps);
209 // compute resulting mean and plot them
210 auto & lr = result.back();
211 if (i < 1 && verbose) Print(lr,"Output for Event" + std::to_string(i));
212 h1->Fill(TMath::Mean(lr.node_data.begin(), lr.node_data.end()));
213 h2->Fill(TMath::Mean(lr.edge_data.begin(), lr.edge_data.end()));
214 h3->Fill(TMath::Mean(lr.global_data.begin(), lr.global_data.end()));
215 }
216 w.Stop();
217 w.Print();
218 check_mem("End evaluation");
219 auto c1 = new TCanvas("c1","SOFIE Results");
220 c1->Divide(1,3);
221 c1->cd(1); h1->Draw();
222 c1->cd(2); h2->Draw();
223 c1->cd(3); h3->Draw();
224
225
226 // compare with the reference
227 auto c2 = new TCanvas("c2","Reference Results");
228 auto file = TFile::Open("graph_data.root");
229 auto o1 = file->Get("h1");
230 auto o2 = file->Get("h2");
231 auto o3 = file->Get("h3");
232 c2->Divide(1,3);
233 c2->cd(1); o1->Draw();
234 c2->cd(2); o2->Draw();
235 c2->cd(3); o3->Draw();
236
237}
#define d(i)
Definition RSha256.hxx:102
#define g(i)
Definition RSha256.hxx:105
#define e(i)
Definition RSha256.hxx:103
#define R__ADD_INCLUDE_PATH(PATH)
Definition Rtypes.h:475
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char filename
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t r
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
R__EXTERN TSystem * gSystem
Definition TSystem.h:582
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
const_iterator begin() const
const_iterator end() const
A "std::vector"-like collection of values implementing handy operation to analyse them.
Definition RVec.hxx:1524
The Canvas class.
Definition TCanvas.h:23
static TFile * Open(const char *name, Option_t *option="", const char *ftitle="", Int_t compress=ROOT::RCompressionSetting::EDefaults::kUseCompiledDefault, Int_t netopt=0)
Create / open a file.
Definition TFile.cxx:3765
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:926
virtual Int_t Fill(Double_t x)
Increment bin with abscissa X by 1.
Definition TH1.cxx:3372
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3076
Random number generator class based on M.
Definition TRandom3.h:27
Stopwatch class.
Definition TStopwatch.h:28
virtual int GetProcInfo(ProcInfo_t *info) const
Returns cpu and memory used by this process into the ProcInfo_t structure.
Definition TSystem.cxx:2500
RVec< Common_t > Concatenate(const RVec< T0 > &v0, const RVec< T1 > &v1)
Return the concatenation of two RVecs.
Definition RVec.hxx:2892
return c1
Definition legend1.C:41
const Int_t n
Definition legend1.C:16
TH1F * h1
Definition legend1.C:5
return c2
Definition legend2.C:14
void Copy(void *source, void *dest)
void Print(std::ostream &os, const OptionType &opt)
std::string ConvertShapeToString(const std::vector< size_t > &shape)
Double_t Mean(Long64_t n, const T *a, const Double_t *w=nullptr)
Returns the weighted mean of an array a with length n.
Definition TMath.h:1176