Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMVA_SOFIE_Inference.py
Go to the documentation of this file.
1### \file
2### \ingroup tutorial_ml
3### \notebook -nodraw
4### This macro provides an example of using a trained model with Keras
5### and make inference using SOFIE directly from Numpy
6### This macro uses as input a Keras model generated with the
7### TMVA_Higgs_Classification.C tutorial
8### You need to run that macro before this one.
9### In this case we are parsing the input file and then run the inference in the same
10### macro making use of the ROOT JITing capability
11###
12###
13### \macro_code
14### \macro_output
15### \author Lorenzo Moneta
16
17import numpy as np
18import ROOT
19
20# check if the input file exists
21modelFile = "Higgs_trained_model.h5"
22if (ROOT.gSystem.AccessPathName(modelFile)) :
23 ROOT.Info("TMVA_SOFIE_RDataFrame","You need to run TMVA_Higgs_Classification.C to generate the Keras trained model")
24 exit()
25
26
27# parse the input Keras model into RModel object
29
30generatedHeaderFile = modelFile.replace(".h5",".hxx")
31print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
32#Generating inference code
34model.OutputGenerated(generatedHeaderFile)
36
37# now compile using ROOT JIT trained model
38modelName = modelFile.replace(".h5","")
39print("compiling SOFIE model ", modelName)
40ROOT.gInterpreter.Declare('#include "' + generatedHeaderFile + '"')
41
42
43generatedHeaderFile = modelFile.replace(".h5",".hxx")
44print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
45#Generating inference
46
47inputFileName = "Higgs_data.root"
48inputFile = str(ROOT.gROOT.GetTutorialDir()) + "/machine_learning/data/" + inputFileName
49
50
51
52
53
54# make SOFIE inference on signal data
55
56df1 = ROOT.RDataFrame("sig_tree", inputFile)
57sigData = df1.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
58#print(sigData)
59
60# stack all the 7 numpy array in a single array (nevents x nvars)
61xsig = np.column_stack(list(sigData.values()))
62dataset_size = xsig.shape[0]
63print("size of data", dataset_size)
64
65#instantiate SOFIE session class
67
68hs = ROOT.TH1D("hs","Signal result",100,0,1)
69for i in range(0,dataset_size):
70 result = session.infer(xsig[i,:])
71 hs.Fill(result[0])
72
73
74# make SOFIE inference on background data
75df2 = ROOT.RDataFrame("bkg_tree", inputFile)
76bkgData = df2.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
77
78xbkg = np.column_stack(list(bkgData.values()))
79dataset_size = xbkg.shape[0]
80
81hb = ROOT.TH1D("hb","Background result",100,0,1)
82for i in range(0,dataset_size):
83 result = session.infer(xbkg[i,:])
84 hb.Fill(result[0])
85
86
87c1 = ROOT.TCanvas()
89hs.SetLineColor("kRed")
90hs.Draw()
91hb.SetLineColor("kBlue")
92hb.Draw("SAME")
94c1.Draw()
95
96
97print("Number of signal entries",hs.GetEntries())
98print("Number of background entries",hb.GetEntries())
99
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...