ROOT
master
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
17
import
numpy
as
np
18
import
ROOT
19
20
# check if the input file exists
21
modelFile =
"Higgs_trained_model.h5"
22
if
(
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
28
model =
ROOT.TMVA.Experimental.SOFIE.PyKeras.Parse
(modelFile)
29
30
generatedHeaderFile =
modelFile.replace
(
".h5"
,
".hxx"
)
31
print(
"Generating inference code for the Keras model from "
,modelFile,
"in the header "
, generatedHeaderFile)
32
#Generating inference code
33
model.Generate
()
34
model.OutputGenerated
(generatedHeaderFile)
35
model.PrintGenerated
()
36
37
# now compile using ROOT JIT trained model
38
modelName =
modelFile.replace
(
".h5"
,
""
)
39
print(
"compiling SOFIE model "
, modelName)
40
ROOT.gInterpreter.Declare
(
'#include "'
+ generatedHeaderFile +
'"'
)
41
42
43
generatedHeaderFile =
modelFile.replace
(
".h5"
,
".hxx"
)
44
print(
"Generating inference code for the Keras model from "
,modelFile,
"in the header "
, generatedHeaderFile)
45
#Generating inference
46
47
inputFileName =
"Higgs_data.root"
48
inputFile = str(
ROOT.gROOT.GetTutorialDir
()) +
"/machine_learning/data/"
+ inputFileName
49
50
51
52
53
54
# make SOFIE inference on signal data
55
56
df1 =
ROOT.RDataFrame
(
"sig_tree"
, inputFile)
57
sigData =
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)
61
xsig =
np.column_stack
(list(
sigData.values
()))
62
dataset_size =
xsig.shape
[0]
63
print(
"size of data"
, dataset_size)
64
65
#instantiate SOFIE session class
66
session =
ROOT.TMVA_SOFIE_Higgs_trained_model.Session
()
67
68
hs =
ROOT.TH1D
(
"hs"
,
"Signal result"
,100,0,1)
69
for
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
75
df2 =
ROOT.RDataFrame
(
"bkg_tree"
, inputFile)
76
bkgData =
df2.AsNumpy
(columns=[
'm_jj'
,
'm_jjj'
,
'm_lv'
,
'm_jlv'
,
'm_bb'
,
'm_wbb'
,
'm_wwbb'
])
77
78
xbkg =
np.column_stack
(list(
bkgData.values
()))
79
dataset_size =
xbkg.shape
[0]
80
81
hb =
ROOT.TH1D
(
"hb"
,
"Background result"
,100,0,1)
82
for
i
in
range
(0,dataset_size):
83
result =
session.infer
(xbkg[i,:])
84
hb.Fill
(result[0])
85
86
87
c1 =
ROOT.TCanvas
()
88
ROOT.gStyle.SetOptStat
(0)
89
hs.SetLineColor
(
"kRed"
)
90
hs.Draw
()
91
hb.SetLineColor
(
"kBlue"
)
92
hb.Draw
(
"SAME"
)
93
c1.BuildLegend
()
94
c1.Draw
()
95
96
97
print(
"Number of signal entries"
,
hs.GetEntries
())
98
print(
"Number of background entries"
,
hb.GetEntries
())
99
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:360
ROOT::Detail::TRangeCast
Definition
TCollection.h:313
ROOT::RDataFrame
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
Definition
RDataFrame.hxx:44
tutorials
machine_learning
TMVA_SOFIE_Inference.py
ROOT master - Reference Guide Generated on Thu Dec 18 2025 14:19:56 (GVA Time) using Doxygen 1.10.0