Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMVA_SOFIE_GNN.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_ml
3## \notebook -nodraw
4##
5## \macro_code
6##
7## \author
8
9import ROOT
10
11# Load system openblas library explicitly if available. This avoids pulling in
12# NumPys builtin openblas later, which will conflict with the system openblas.
13ROOT.gInterpreter.Load("libopenblaso.so")
14
15import time
16
17import graph_nets as gn
18import numpy as np
19import sonnet as snt
20from graph_nets import utils_tf
21
22# defining graph properties
23num_nodes = 5
24num_edges = 20
25snd = np.array([1, 2, 3, 4, 2, 3, 4, 3, 4, 4, 0, 0, 0, 0, 1, 1, 1, 2, 2, 3], dtype="int32")
26rec = np.array([0, 0, 0, 0, 1, 1, 1, 2, 2, 3, 1, 2, 3, 4, 2, 3, 4, 3, 4, 4], dtype="int32")
27node_size = 4
28edge_size = 4
29global_size = 1
30LATENT_SIZE = 100
31NUM_LAYERS = 4
32processing_steps = 5
33
34
35# method for returning dictionary of graph data
36def get_graph_data_dict(num_nodes, num_edges, NODE_FEATURE_SIZE=2, EDGE_FEATURE_SIZE=2, GLOBAL_FEATURE_SIZE=1):
37 return {
38 "globals": 10 * np.random.rand(GLOBAL_FEATURE_SIZE).astype(np.float32) - 5.0,
39 "nodes": 10 * np.random.rand(num_nodes, NODE_FEATURE_SIZE).astype(np.float32) - 5.0,
40 "edges": 10 * np.random.rand(num_edges, EDGE_FEATURE_SIZE).astype(np.float32) - 5.0,
41 "senders": snd,
42 "receivers": rec,
43 }
44
45
46# method to instantiate mlp model to be added in GNN
47def make_mlp_model():
48 return snt.Sequential(
49 [
50 snt.nets.MLP([LATENT_SIZE] * NUM_LAYERS, activate_final=True),
51 snt.LayerNorm(axis=-1, create_offset=True, create_scale=True),
52 ]
53 )
54
55
56# defining GraphIndependent class with MLP edge, node, and global models.
58 def __init__(self, name="MLPGraphIndependent"):
59 super(MLPGraphIndependent, self).__init__(name=name)
60 self._network = gn.modules.GraphIndependent(
61 edge_model_fn=lambda: snt.nets.MLP([LATENT_SIZE] * NUM_LAYERS, activate_final=True),
62 node_model_fn=lambda: snt.nets.MLP([LATENT_SIZE] * NUM_LAYERS, activate_final=True),
63 global_model_fn=lambda: snt.nets.MLP([LATENT_SIZE] * NUM_LAYERS, activate_final=True),
64 )
65
66 def __call__(self, inputs):
67 return self._network(inputs)
68
69
70# defining Graph network class with MLP edge, node, and global models.
72 def __init__(self, name="MLPGraphNetwork"):
73 super(MLPGraphNetwork, self).__init__(name=name)
74 self._network = gn.modules.GraphNetwork(
75 edge_model_fn=make_mlp_model, node_model_fn=make_mlp_model, global_model_fn=make_mlp_model
76 )
77
78 def __call__(self, inputs):
79 return self._network(inputs)
80
81
82# defining a Encode-Process-Decode module for LHCb toy model
84
85 def __init__(self, name="EncodeProcessDecode"):
86 super(EncodeProcessDecode, self).__init__(name=name)
87 self._encoder = MLPGraphIndependent()
88 self._core = MLPGraphNetwork()
89 self._decoder = MLPGraphIndependent()
90 self._output_transform = MLPGraphIndependent()
91
92 def __call__(self, input_op, num_processing_steps):
93 latent = self._encoder(input_op)
94 latent0 = latent
95 output_ops = []
96 for _ in range(num_processing_steps):
97 core_input = utils_tf.concat([latent0, latent], axis=1)
98 latent = self._core(core_input)
99 decoded_op = self._decoder(latent)
100 output_ops.append(self._output_transform(decoded_op))
101 return output_ops
102
103
104# Instantiating EncodeProcessDecode Model
105ep_model = EncodeProcessDecode()
106
107# Initializing randomized input data
108GraphData = get_graph_data_dict(num_nodes, num_edges, node_size, edge_size, global_size)
109
110# input_graphs is a tuple representing the initial data
111input_graph_data = utils_tf.data_dicts_to_graphs_tuple([GraphData])
112
113# Initializing randomized input data for core
114# note that the core network has as input a double number of features
115CoreGraphData = get_graph_data_dict(num_nodes, num_edges, 2 * LATENT_SIZE, 2 * LATENT_SIZE, 2 * LATENT_SIZE)
116input_core_graph_data = utils_tf.data_dicts_to_graphs_tuple([CoreGraphData])
117
118# initialize graph data for decoder (input is LATENT_SIZE)
119DecodeGraphData = get_graph_data_dict(num_nodes, num_edges, LATENT_SIZE, LATENT_SIZE, LATENT_SIZE)
120
121# Make prediction of GNN
122output_gn = ep_model(input_graph_data, processing_steps)
123# print("---> Input:\n",input_graph_data)
124# print("\n\n------> Input core data:\n",input_core_graph_data)
125# print("\n\n---> Output:\n",output_gn)
126
127# Make SOFIE Model
129 ep_model._encoder._network, GraphData, filename="gnn_encoder"
130)
133
135 ep_model._core._network, CoreGraphData, filename="gnn_core"
136)
139
141 ep_model._decoder._network, DecodeGraphData, filename="gnn_decoder"
142)
145
147 ep_model._output_transform._network, DecodeGraphData, filename="gnn_output_transform"
148)
151
152# Compile now the generated C++ code from SOFIE
153gen_code = '''#pragma cling optimize(2)
154#include "gnn_encoder.hxx"
155#include "gnn_core.hxx"
156#include "gnn_decoder.hxx"
157#include "gnn_output_transform.hxx"'''
159
160
161# helper function to print SOFIE GNN data structure
162def PrintSofie(output, printShape=False):
166 if printShape:
167 print("SOFIE data ... shapes", n.shape, e.shape, g.shape)
168 print(
169 " node data",
170 n.reshape(
171 n.size,
172 ),
173 )
174 print(
175 " edge data",
176 e.reshape(
177 e.size,
178 ),
179 )
180 print(
181 " global data",
182 g.reshape(
183 g.size,
184 ),
185 )
186
187
188def CopyData(input_data):
189 output_data = ROOT.TMVA.Experimental.SOFIE.Copy(input_data)
190 return output_data
191
192
193# Build SOFIE GNN Model and run inference
194class SofieGNN:
195 def __init__(self):
196 self.encoder_session = ROOT.TMVA_SOFIE_gnn_encoder.Session()
197 self.core_session = ROOT.TMVA_SOFIE_gnn_core.Session()
198 self.decoder_session = ROOT.TMVA_SOFIE_gnn_decoder.Session()
199 self.output_transform_session = ROOT.TMVA_SOFIE_gnn_output_transform.Session()
200
201 def infer(self, graphData):
202 # copy the input data
203 input_data = CopyData(graphData)
204
205 # running inference on sofie
206 self.encoder_session.infer(input_data)
207 latent0 = CopyData(input_data)
208 latent = input_data
209 output_ops = []
210 for _ in range(processing_steps):
211 core_input = ROOT.TMVA.Experimental.SOFIE.Concatenate(latent0, latent, axis=1)
212 self.core_session.infer(core_input)
213 latent = CopyData(core_input)
214 self.decoder_session.infer(core_input)
215 self.output_transform_session.infer(core_input)
216 output = CopyData(core_input)
217 output_ops.append(output)
218
219 return output_ops
220
221
222# Test both GNN on some simulated events
223def GenerateData():
224 data = get_graph_data_dict(num_nodes, num_edges, node_size, edge_size, global_size)
225 return data
226
227
228numevts = 40
229dataSet = []
230for i in range(0, numevts):
231 data = GenerateData()
232 dataSet.append(data)
233
234# Run graph_nets model
235# First we convert input data to the required input format
236gnetData = []
237for i in range(0, numevts):
238 graphData = dataSet[i]
239 gnet_data_i = utils_tf.data_dicts_to_graphs_tuple([graphData])
240 gnetData.append(gnet_data_i)
241
242
243# Function to run the graph net
244def RunGNet(inputGraphData):
245 output_gn = ep_model(inputGraphData, processing_steps)
246 return output_gn
247
248
249start = time.time()
250hG = ROOT.TH1D("hG", "Result from graphnet", 20, 1, 0)
251for i in range(0, numevts):
252 out = RunGNet(gnetData[i])
253 g = out[1].globals.numpy()
254 hG.Fill(np.mean(g))
255
256end = time.time()
257print("elapsed time for ", numevts, "events = ", end - start)
258
259# running SOFIE-GNN
260sofieData = []
261for i in range(0, numevts):
262 graphData = dataSet[i]
267 input_data.edge_index = ROOT.TMVA.Experimental.AsRTensor(np.stack((graphData["receivers"], graphData["senders"])))
268 sofieData.append(input_data)
269
270
271endSC = time.time()
272print("time to convert data to SOFIE format", endSC - end)
273
274hS = ROOT.TH1D("hS", "Result from SOFIE", 20, 1, 0)
275start0 = time.time()
276gnn = SofieGNN()
277start = time.time()
278print("time to create SOFIE GNN class", start - start0)
279for i in range(0, numevts):
280 # print("inference event....",i)
281 out = gnn.infer(sofieData[i])
282 g = np.asarray(out[1].global_data)
283 hS.Fill(np.mean(g))
284
285end = time.time()
286print("elapsed time for ", numevts, "events = ", end - start)
287
288c0 = ROOT.TCanvas()
289c0.Divide(1, 2)
290c1 = c0.cd(1)
291c1.Divide(2, 1)
292c1.cd(1)
293hG.Draw()
294c1.cd(2)
295hS.Draw()
296
297hDe = ROOT.TH1D("hDe", "Difference for edge data", 40, 1, 0)
298hDn = ROOT.TH1D("hDn", "Difference for node data", 40, 1, 0)
299hDg = ROOT.TH1D("hDg", "Difference for global data", 40, 1, 0)
300# compute differences between SOFIE and GNN
301for i in range(0, numevts):
302 outSofie = gnn.infer(sofieData[i])
303 outGnet = RunGNet(gnetData[i])
304 edgesG = outGnet[1].edges.numpy()
305 edgesS = np.asarray(outSofie[1].edge_data)
306 if i == 0:
307 print(edgesG.shape)
308 for j in range(0, edgesG.shape[0]):
309 for k in range(0, edgesG.shape[1]):
310 hDe.Fill(edgesG[j, k] - edgesS[j, k])
311
312 nodesG = outGnet[1].nodes.numpy()
313 nodesS = np.asarray(outSofie[1].node_data)
314 for j in range(0, nodesG.shape[0]):
315 for k in range(0, nodesG.shape[1]):
316 hDn.Fill(nodesG[j, k] - nodesS[j, k])
317
318 globG = outGnet[1].globals.numpy()
319 globS = np.asarray(outSofie[1].global_data)
320 for j in range(0, globG.shape[1]):
321 hDg.Fill(globG[0, j] - globS[j])
322
323
324c2 = c0.cd(2)
325c2.Divide(3, 1)
326c2.cd(1)
327hDe.Draw()
328c2.cd(2)
329hDn.Draw()
330c2.cd(3)
331hDg.Draw()
332
333c0.Draw()
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.