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