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
17from os.path import exists
18
19import numpy as np
20import ROOT
21
22# check if the input file exists
23modelFile = "HiggsModel.keras"
24
25if not exists(modelFile):
26 raise FileNotFoundError("You need to run TMVA_Higgs_Classification.C to generate the Keras trained model")
27
28
29# parse the input Keras model into RModel object
31
32generatedHeaderFile = modelFile.replace(".keras",".hxx")
33print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
34#Generating inference code
36model.OutputGenerated(generatedHeaderFile)
38
39# now compile using ROOT JIT trained model
40modelName = modelFile.replace(".keras","")
41print("compiling SOFIE model ", modelName)
42ROOT.gInterpreter.Declare('#include "' + generatedHeaderFile + '"')
43
44
45generatedHeaderFile = modelFile.replace(".keras",".hxx")
46print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
47#Generating inference
48
49inputFileName = "Higgs_data.root"
50inputFile = str(ROOT.gROOT.GetTutorialDir()) + "/machine_learning/data/" + inputFileName
51
52
53
54
55
56# make SOFIE inference on signal data
57
58df1 = ROOT.RDataFrame("sig_tree", inputFile)
59sigData = df1.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
60#print(sigData)
61
62# stack all the 7 numpy array in a single array (nevents x nvars)
63xsig = np.column_stack(list(sigData.values()))
64dataset_size = xsig.shape[0]
65print("size of signal data", dataset_size)
66
67#instantiate SOFIE session class
68#session = ROOT.TMVA_SOFIE_HiggsModel.Session()
69#get the sofie session namespace
70sofie = getattr(ROOT, 'TMVA_SOFIE_' + modelName)
71session = sofie.Session()
72
73print("Evaluating SOFIE models on signal data")
74hs = ROOT.TH1D("hs","Signal result",100,0,1)
75for i in range(0,dataset_size):
76 result = session.infer(xsig[i,:])
77 if (i % dataset_size/10 == 0) :
78 print("result for signal event ",i,result[0])
79 hs.Fill(result[0])
80
81print("using RDsataFrame to extract input data in a numpy array")
82# make SOFIE inference on background data
83df2 = ROOT.RDataFrame("bkg_tree", inputFile)
84bkgData = df2.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
85
86xbkg = np.column_stack(list(bkgData.values()))
87dataset_size = xbkg.shape[0]
88print("size of background data", dataset_size)
89
90hb = ROOT.TH1D("hb","Background result",100,0,1)
91for i in range(0,dataset_size):
92 result = session.infer(xbkg[i,:])
93 if (i % dataset_size/10 == 0) :
94 print("result for background event ",i,result[0])
95
96 hb.Fill(result[0])
97
98
99c1 = ROOT.TCanvas()
101hs.SetLineColor("kRed")
102hs.Draw()
103hb.SetLineColor("kBlue")
104hb.Draw("SAME")
106c1.Draw()
107
108
109print("Number of signal entries",hs.GetEntries())
110print("Number of background entries",hb.GetEntries())
111
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 ,...