Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMVA_SOFIE_GNN.py
Go to the documentation of this file.
1import ROOT
2
3import numpy as np
4import graph_nets as gn
5from graph_nets import utils_tf
6import sonnet as snt
7import time
8
9# defining graph properties
10num_nodes=5
11num_edges=20
12snd = np.array([1,2,3,4,2,3,4,3,4,4,0,0,0,0,1,1,1,2,2,3], dtype='int32')
13rec = np.array([0,0,0,0,1,1,1,2,2,3,1,2,3,4,2,3,4,3,4,4], dtype='int32')
14node_size=4
15edge_size=4
16global_size=1
17LATENT_SIZE = 100
18NUM_LAYERS = 4
19processing_steps = 5
20
21# method for returning dictionary of graph data
22def get_graph_data_dict(num_nodes, num_edges, NODE_FEATURE_SIZE=2, EDGE_FEATURE_SIZE=2, GLOBAL_FEATURE_SIZE=1):
23 return {
24 "globals": 10*np.random.rand(GLOBAL_FEATURE_SIZE).astype(np.float32)-5.,
25 "nodes": 10*np.random.rand(num_nodes, NODE_FEATURE_SIZE).astype(np.float32)-5.,
26 "edges": 10*np.random.rand(num_edges, EDGE_FEATURE_SIZE).astype(np.float32)-5.,
27 "senders": snd,
28 "receivers": rec
29 }
30
31# method to instantiate mlp model to be added in GNN
33 return snt.Sequential([
34 snt.nets.MLP([LATENT_SIZE]*NUM_LAYERS, activate_final=True),
35 snt.LayerNorm(axis=-1, create_offset=True, create_scale=True)
36 ])
37
38# defining GraphIndependent class with MLP edge, node, and global models.
39class MLPGraphIndependent(snt.Module):
40 def __init__(self, name="MLPGraphIndependent"):
41 super(MLPGraphIndependent, self).__init__(name=name)
42 self._network = gn.modules.GraphIndependent(
43 edge_model_fn = lambda: snt.nets.MLP([LATENT_SIZE]*NUM_LAYERS, activate_final=True),
44 node_model_fn = lambda: snt.nets.MLP([LATENT_SIZE]*NUM_LAYERS, activate_final=True),
45 global_model_fn = lambda: snt.nets.MLP([LATENT_SIZE]*NUM_LAYERS, activate_final=True))
46
47 def __call__(self, inputs):
48 return self._network(inputs)
49
50# defining Graph network class with MLP edge, node, and global models.
51class MLPGraphNetwork(snt.Module):
52 def __init__(self, name="MLPGraphNetwork"):
53 super(MLPGraphNetwork, self).__init__(name=name)
54 self._network = gn.modules.GraphNetwork(
55 edge_model_fn=make_mlp_model,
56 node_model_fn=make_mlp_model,
57 global_model_fn=make_mlp_model)
58
59 def __call__(self, inputs):
60 return self._network(inputs)
61
62# defining a Encode-Process-Decode module for LHCb toy model
63class EncodeProcessDecode(snt.Module):
64
65 def __init__(self,
66 name="EncodeProcessDecode"):
67 super(EncodeProcessDecode, self).__init__(name=name)
72
73 def __call__(self, input_op, num_processing_steps):
74 latent = self._encoder(input_op)
75 latent0 = latent
76 output_ops = []
77 for _ in range(num_processing_steps):
78 core_input = utils_tf.concat([latent0, latent], axis=1)
79 latent = self._core(core_input)
80 decoded_op = self._decoder(latent)
81 output_ops.append(self._output_transform(decoded_op))
82 return output_ops
83
84
85# Instantiating EncodeProcessDecode Model
87
88# Initializing randomized input data
89GraphData = get_graph_data_dict(num_nodes,num_edges, node_size, edge_size, global_size)
90
91#input_graphs is a tuple representing the initial data
92input_graph_data = utils_tf.data_dicts_to_graphs_tuple([GraphData])
93
94# Initializing randomized input data for core
95# note that the core network has as input a double number of features
96CoreGraphData = get_graph_data_dict(num_nodes, num_edges, 2*LATENT_SIZE, 2*LATENT_SIZE, 2*LATENT_SIZE)
97input_core_graph_data = utils_tf.data_dicts_to_graphs_tuple([CoreGraphData])
98
99#initialize graph data for decoder (input is LATENT_SIZE)
100DecodeGraphData = get_graph_data_dict(num_nodes,num_edges, LATENT_SIZE, LATENT_SIZE, LATENT_SIZE)
101
102# Make prediction of GNN
103output_gn = ep_model(input_graph_data, processing_steps)
104print("---> Input:\n",input_graph_data)
105print("\n\n------> Input core data:\n",input_core_graph_data)
106print("\n\n---> Output:\n",output_gn)
107
108# Make SOFIE Model
109encoder = ROOT.TMVA.Experimental.SOFIE.RModel_GraphIndependent.ParseFromMemory(ep_model._encoder._network, GraphData, filename = "encoder")
110encoder.Generate()
111encoder.OutputGenerated()
112
113core = ROOT.TMVA.Experimental.SOFIE.RModel_GNN.ParseFromMemory(ep_model._core._network, CoreGraphData, filename = "core")
114core.Generate()
115core.OutputGenerated()
116
117decoder = ROOT.TMVA.Experimental.SOFIE.RModel_GraphIndependent.ParseFromMemory(ep_model._decoder._network, DecodeGraphData, filename = "decoder")
118decoder.Generate()
119decoder.OutputGenerated()
120
121output_transform = ROOT.TMVA.Experimental.SOFIE.RModel_GraphIndependent.ParseFromMemory(ep_model._output_transform._network, DecodeGraphData, filename = "output_transform")
122output_transform.Generate()
123output_transform.OutputGenerated()
124
125# Compile now the generated C++ code from SOFIE
126ROOT.gInterpreter.Declare('#pragma cling optimize(2)')
127ROOT.gInterpreter.Declare('#include "encoder.hxx"')
128ROOT.gInterpreter.Declare('#include "core.hxx"')
129ROOT.gInterpreter.Declare('#include "decoder.hxx"')
130ROOT.gInterpreter.Declare('#include "output_transform.hxx"')
131
132#helper function to print SOFIE GNN data structure
133def PrintSofie(output, printShape = False):
134 n = np.asarray(output.node_data)
135 e = np.asarray(output.edge_data)
136 g = np.asarray(output.global_data)
137 if (printShape) :
138 print("SOFIE data ... shapes",n.shape,e.shape,g.shape)
139 print(" node data", n.reshape(n.size,))
140 print(" edge data", e.reshape(e.size,))
141 print(" global data",g.reshape(g.size,))
142
143def CopyData(input_data) :
144 output_data = ROOT.TMVA.Experimental.SOFIE.Copy(input_data)
145 return output_data
146
147# Build SOFIE GNN Model and run inference
149 def __init__(self):
150 self.encoder_session = ROOT.TMVA_SOFIE_encoder.Session()
151 self.core_session = ROOT.TMVA_SOFIE_core.Session()
152 self.decoder_session = ROOT.TMVA_SOFIE_decoder.Session()
153 self.output_transform_session = ROOT.TMVA_SOFIE_output_transform.Session()
154
155 def infer(self, graphData):
156 # copy the input data
157 input_data = CopyData(graphData)
158
159 # running inference on sofie
160 self.encoder_session.infer(input_data)
161 latent0 = CopyData(input_data)
162 latent = input_data
163 output_ops = []
164 for _ in range(processing_steps):
165 core_input = ROOT.TMVA.Experimental.SOFIE.Concatenate(latent0, latent, axis=1)
166 self.core_session.infer(core_input)
167 latent = CopyData(core_input)
168 self.decoder_session.infer(core_input)
169 self.output_transform_session.infer(core_input)
170 output = CopyData(core_input)
171 output_ops.append(output)
172
173 return output_ops
174
175# Test both GNN on some simulated events
177 data = get_graph_data_dict(num_nodes,num_edges, node_size, edge_size, global_size)
178 return data
179
180numevts = 100
181dataSet = []
182for i in range(0,numevts):
184 dataSet.append(data)
185
186# Run graph_nets model
187# First we convert input data to the required input format
188gnetData = []
189for i in range(0,numevts):
190 graphData = dataSet[i]
191 gnet_data_i = utils_tf.data_dicts_to_graphs_tuple([graphData])
192 gnetData.append(gnet_data_i)
193
194# Function to run the graph net
195def RunGNet(inputGraphData) :
196 output_gn = ep_model(inputGraphData, processing_steps)
197 return output_gn
198
199start = time.time()
200hG = ROOT.TH1D("hG","Result from graphnet",100,1,0)
201for i in range(0,numevts):
202 out = RunGNet(gnetData[i])
203 g = out[1].globals.numpy()
204 hG.Fill(np.mean(g))
205
206end = time.time()
207print("elapsed time for ",numevts,"events = ",end-start)
208
209# running SOFIE-GNN
210sofieData = []
211for i in range(0,numevts):
212 graphData = dataSet[i]
213 input_data = ROOT.TMVA.Experimental.SOFIE.GNN_Data()
214 input_data.node_data = ROOT.TMVA.Experimental.AsRTensor(graphData['nodes'])
215 input_data.edge_data = ROOT.TMVA.Experimental.AsRTensor(graphData['edges'])
216 input_data.global_data = ROOT.TMVA.Experimental.AsRTensor(graphData['globals'])
217 #make sure dtype of graphData['receivers'] and senders is int32
218 input_data.receivers = graphData['receivers']
219 input_data.senders = graphData['senders']
220 sofieData.append(input_data)
221
222print("SOFIE Data: first event")
223print("receivers",sofieData[0].receivers)
224print("senders",sofieData[0].senders)
225
226endSC = time.time()
227print("time to convert data to SOFIE format",endSC-end)
228
229hS = ROOT.TH1D("hS","Result from SOFIE",100,1,0)
230start0 = time.time()
231gnn = SofieGNN()
232start = time.time()
233print("time to create SOFIE GNN class", start-start0)
234for i in range(0,numevts):
235 #print("inference event....",i)
236 out = gnn.infer(sofieData[i])
237 g = np.asarray(out[1].global_data)
238 hS.Fill(np.mean(g))
239
240end = time.time()
241print("elapsed time for ",numevts,"events = ",end-start)
242
243c0 = ROOT.TCanvas()
244c0.Divide(1,2)
245c1 = c0.cd(1)
246c1.Divide(2,1)
247c1.cd(1)
248hG.Draw()
249c1.cd(2)
250hS.Draw()
251
252hDe = ROOT.TH1D("hDe","Difference for edge data",100,1,0)
253hDn = ROOT.TH1D("hDn","Difference for node data",100,1,0)
254hDg = ROOT.TH1D("hDg","Difference for global data",100,1,0)
255#compute differences between SOFIE and GNN
256for i in range(0,numevts):
257 outSofie = gnn.infer(sofieData[i])
258 outGnet = RunGNet(gnetData[i])
259 edgesG = outGnet[1].edges.numpy()
260 edgesS = np.asarray(outSofie[1].edge_data)
261 if (i == 0) : print(edgesG.shape)
262 for j in range(0,edgesG.shape[0]) :
263 for k in range(0,edgesG.shape[1]) :
264 hDe.Fill(edgesG[j,k]-edgesS[j,k])
265
266 nodesG = outGnet[1].nodes.numpy()
267 nodesS = np.asarray(outSofie[1].node_data)
268 for j in range(0,nodesG.shape[0]) :
269 for k in range(0,nodesG.shape[1]) :
270 hDn.Fill(nodesG[j,k]-nodesS[j,k])
271
272 globG = outGnet[1].globals.numpy()
273 globS = np.asarray(outSofie[1].global_data)
274 for j in range(0,globG.shape[1]) :
275 hDg.Fill(globG[0,j]-globS[j])
276
277
278c2 = c0.cd(2)
279c2.Divide(3,1)
280c2.cd(1)
281hDe.Draw()
282c2.cd(2)
283hDn.Draw()
284c2.cd(3)
285hDg.Draw()
286
287c0.Draw()
288
__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)