Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
TMVA_SOFIE_Inference.py File Reference

Detailed Description

View in nbviewer Open in SWAN
This macro provides an example of using a trained model with Keras and make inference using SOFIE directly from Numpy This macro uses as input a Keras model generated with the TMVA_Higgs_Classification.C tutorial You need to run that macro before this one.

In this case we are parsing the input file and then run the inference in the same macro making use of the ROOT JITing capability

import ROOT
import numpy as np
ROOT.TMVA.PyMethodBase.PyInitialize()
# check if the input file exists
modelFile = "Higgs_trained_model.h5"
if (ROOT.gSystem.AccessPathName(modelFile)) :
ROOT.Info("TMVA_SOFIE_RDataFrame","You need to run TMVA_Higgs_Classification.C to generate the Keras trained model")
exit()
# parse the input Keras model into RModel object
model = ROOT.TMVA.Experimental.SOFIE.PyKeras.Parse(modelFile)
generatedHeaderFile = modelFile.replace(".h5",".hxx")
print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
#Generating inference code
model.Generate()
model.OutputGenerated(generatedHeaderFile)
model.PrintGenerated()
# now compile using ROOT JIT trained model
modelName = modelFile.replace(".h5","")
print("compiling SOFIE model ", modelName)
ROOT.gInterpreter.Declare('#include "' + generatedHeaderFile + '"')
generatedHeaderFile = modelFile.replace(".h5",".hxx")
print("Generating inference code for the Keras model from ",modelFile,"in the header ", generatedHeaderFile)
#Generating inference
inputFileName = "Higgs_data.root"
inputFile = "http://root.cern.ch/files/" + inputFileName
# make SOFIE inference on signal data
df1 = ROOT.RDataFrame("sig_tree", inputFile)
sigData = df1.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
#print(sigData)
# stack all the 7 numpy array in a single array (nevents x nvars)
xsig = np.column_stack(list(sigData.values()))
dataset_size = xsig.shape[0]
print("size of data", dataset_size)
#instantiate SOFIE session class
session = ROOT.TMVA_SOFIE_Higgs_trained_model.Session()
hs = ROOT.TH1D("hs","Signal result",100,0,1)
for i in range(0,dataset_size):
result = session.infer(xsig[i,:])
hs.Fill(result[0])
# make SOFIE inference on background data
df2 = ROOT.RDataFrame("bkg_tree", inputFile)
bkgData = df2.AsNumpy(columns=['m_jj', 'm_jjj', 'm_lv', 'm_jlv', 'm_bb', 'm_wbb', 'm_wwbb'])
xbkg = np.column_stack(list(bkgData.values()))
dataset_size = xbkg.shape[0]
hb = ROOT.TH1D("hb","Background result",100,0,1)
for i in range(0,dataset_size):
result = session.infer(xbkg[i,:])
hb.Fill(result[0])
c1 = ROOT.TCanvas()
ROOT.gStyle.SetOptStat(0)
hs.SetLineColor(ROOT.kRed)
hs.Draw()
hb.SetLineColor(ROOT.kBlue)
hb.Draw("SAME")
c1.BuildLegend()
c1.Draw()
print("Number of signal entries",hs.GetEntries())
print("Number of background entries",hb.GetEntries())
ROOT's RDataFrame offers a modern, high-level interface for analysis of data stored in TTree ,...
Model has not a defined batch size assume is 1 - input shape for tensor dense_input : { 1 , 7 }
//Code generated automatically by TMVA for Inference of Model file [Higgs_trained_model.h5] at [Thu Dec 19 08:55:19 2024]
#ifndef ROOT_TMVA_SOFIE_HIGGS_TRAINED_MODEL
#define ROOT_TMVA_SOFIE_HIGGS_TRAINED_MODEL
#include <algorithm>
#include <vector>
#include <cmath>
#include "TMVA/SOFIE_common.hxx"
#include <fstream>
namespace TMVA_SOFIE_Higgs_trained_model{
namespace BLAS{
extern "C" void sgemv_(const char * trans, const int * m, const int * n, const float * alpha, const float * A,
const int * lda, const float * X, const int * incx, const float * beta, const float * Y, const int * incy);
extern "C" void sgemm_(const char * transa, const char * transb, const int * m, const int * n, const int * k,
const float * alpha, const float * A, const int * lda, const float * B, const int * ldb,
const float * beta, float * C, const int * ldc);
}//BLAS
struct Session {
// initialized tensors
std::vector<float> fTensor_dense_4kernel0 = std::vector<float>(128);
float * tensor_dense_4kernel0 = fTensor_dense_4kernel0.data();
std::vector<float> fTensor_dense_3bias0 = std::vector<float>(64);
float * tensor_dense_3bias0 = fTensor_dense_3bias0.data();
std::vector<float> fTensor_dense_3kernel0 = std::vector<float>(4096);
float * tensor_dense_3kernel0 = fTensor_dense_3kernel0.data();
std::vector<float> fTensor_dense_2kernel0 = std::vector<float>(4096);
float * tensor_dense_2kernel0 = fTensor_dense_2kernel0.data();
std::vector<float> fTensor_dense_1bias0 = std::vector<float>(64);
float * tensor_dense_1bias0 = fTensor_dense_1bias0.data();
std::vector<float> fTensor_dense_4bias0 = std::vector<float>(2);
float * tensor_dense_4bias0 = fTensor_dense_4bias0.data();
std::vector<float> fTensor_dense_1kernel0 = std::vector<float>(4096);
float * tensor_dense_1kernel0 = fTensor_dense_1kernel0.data();
std::vector<float> fTensor_densebias0 = std::vector<float>(64);
float * tensor_densebias0 = fTensor_densebias0.data();
std::vector<float> fTensor_dense_2bias0 = std::vector<float>(64);
float * tensor_dense_2bias0 = fTensor_dense_2bias0.data();
std::vector<float> fTensor_densekernel0 = std::vector<float>(448);
float * tensor_densekernel0 = fTensor_densekernel0.data();
//--- declare and allocate the intermediate tensors
std::vector<float> fTensor_dense_4Sigmoid0 = std::vector<float>(2);
float * tensor_dense_4Sigmoid0 = fTensor_dense_4Sigmoid0.data();
std::vector<float> fTensor_dense_4Dense = std::vector<float>(2);
float * tensor_dense_4Dense = fTensor_dense_4Dense.data();
std::vector<float> fTensor_densebias0bcast = std::vector<float>(64);
float * tensor_densebias0bcast = fTensor_densebias0bcast.data();
std::vector<float> fTensor_dense_1bias0bcast = std::vector<float>(64);
float * tensor_dense_1bias0bcast = fTensor_dense_1bias0bcast.data();
std::vector<float> fTensor_dense_3bias0bcast = std::vector<float>(64);
float * tensor_dense_3bias0bcast = fTensor_dense_3bias0bcast.data();
std::vector<float> fTensor_dense_1Dense = std::vector<float>(64);
float * tensor_dense_1Dense = fTensor_dense_1Dense.data();
std::vector<float> fTensor_dense_1Relu0 = std::vector<float>(64);
float * tensor_dense_1Relu0 = fTensor_dense_1Relu0.data();
std::vector<float> fTensor_dense_2Dense = std::vector<float>(64);
float * tensor_dense_2Dense = fTensor_dense_2Dense.data();
std::vector<float> fTensor_denseRelu0 = std::vector<float>(64);
float * tensor_denseRelu0 = fTensor_denseRelu0.data();
std::vector<float> fTensor_dense_2Relu0 = std::vector<float>(64);
float * tensor_dense_2Relu0 = fTensor_dense_2Relu0.data();
std::vector<float> fTensor_dense_3Dense = std::vector<float>(64);
float * tensor_dense_3Dense = fTensor_dense_3Dense.data();
std::vector<float> fTensor_denseDense = std::vector<float>(64);
float * tensor_denseDense = fTensor_denseDense.data();
std::vector<float> fTensor_dense_3Relu0 = std::vector<float>(64);
float * tensor_dense_3Relu0 = fTensor_dense_3Relu0.data();
std::vector<float> fTensor_dense_2bias0bcast = std::vector<float>(64);
float * tensor_dense_2bias0bcast = fTensor_dense_2bias0bcast.data();
std::vector<float> fTensor_dense_4bias0bcast = std::vector<float>(2);
float * tensor_dense_4bias0bcast = fTensor_dense_4bias0bcast.data();
Session(std::string filename ="Higgs_trained_model.dat") {
//--- reading weights from file
std::ifstream f;
f.open(filename);
if (!f.is_open()) {
throw std::runtime_error("tmva-sofie failed to open file " + filename + " for input weights");
}
std::string tensor_name;
size_t length;
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_4kernel0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_4kernel0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 128) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 128 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_4kernel0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_4kernel0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_3bias0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_3bias0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 64) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 64 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_3bias0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_3bias0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_3kernel0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_3kernel0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 4096) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 4096 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_3kernel0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_3kernel0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_2kernel0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_2kernel0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 4096) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 4096 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_2kernel0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_2kernel0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_1bias0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_1bias0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 64) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 64 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_1bias0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_1bias0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_4bias0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_4bias0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 2) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 2 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_4bias0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_4bias0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_1kernel0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_1kernel0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 4096) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 4096 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_1kernel0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_1kernel0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_densebias0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_densebias0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 64) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 64 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_densebias0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_densebias0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_dense_2bias0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_dense_2bias0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 64) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 64 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_dense_2bias0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_dense_2bias0");
}
f >> tensor_name >> length;
if (tensor_name != "tensor_densekernel0" ) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor name; expected name is tensor_densekernel0 , read " + tensor_name;
throw std::runtime_error(err_msg);
}
if (length != 448) {
std::string err_msg = "TMVA-SOFIE failed to read the correct tensor size; expected size is 448 , read " + std::to_string(length) ;
throw std::runtime_error(err_msg);
}
for (size_t i = 0; i < length; ++i)
f >> tensor_densekernel0[i];
if (f.fail()) {
throw std::runtime_error("TMVA-SOFIE failed to read the values for tensor tensor_densekernel0");
}
f.close();
//---- allocate the intermediate dynamic tensors
//--- broadcast bias tensor densebias0for Gemm op
{
float * data = TMVA::Experimental::SOFIE::UTILITY::UnidirectionalBroadcast<float>(tensor_densebias0,{ 64 }, { 1 , 64 });
std::copy(data, data + 64, tensor_densebias0bcast);
delete [] data;
}
//--- broadcast bias tensor dense_1bias0for Gemm op
{
float * data = TMVA::Experimental::SOFIE::UTILITY::UnidirectionalBroadcast<float>(tensor_dense_1bias0,{ 64 }, { 1 , 64 });
std::copy(data, data + 64, tensor_dense_1bias0bcast);
delete [] data;
}
//--- broadcast bias tensor dense_2bias0for Gemm op
{
float * data = TMVA::Experimental::SOFIE::UTILITY::UnidirectionalBroadcast<float>(tensor_dense_2bias0,{ 64 }, { 1 , 64 });
std::copy(data, data + 64, tensor_dense_2bias0bcast);
delete [] data;
}
//--- broadcast bias tensor dense_3bias0for Gemm op
{
float * data = TMVA::Experimental::SOFIE::UTILITY::UnidirectionalBroadcast<float>(tensor_dense_3bias0,{ 64 }, { 1 , 64 });
std::copy(data, data + 64, tensor_dense_3bias0bcast);
delete [] data;
}
//--- broadcast bias tensor dense_4bias0for Gemm op
{
float * data = TMVA::Experimental::SOFIE::UTILITY::UnidirectionalBroadcast<float>(tensor_dense_4bias0,{ 2 }, { 1 , 2 });
std::copy(data, data + 2, tensor_dense_4bias0bcast);
delete [] data;
}
}
std::vector<float> infer(float* tensor_dense_input){
//--------- Gemm
char op_0_transA = 'n';
char op_0_transB = 'n';
int op_0_m = 1;
int op_0_n = 64;
int op_0_k = 7;
float op_0_alpha = 1;
float op_0_beta = 1;
int op_0_lda = 7;
int op_0_ldb = 64;
std::copy(tensor_densebias0bcast, tensor_densebias0bcast + 64, tensor_denseDense);
BLAS::sgemm_(&op_0_transB, &op_0_transA, &op_0_n, &op_0_m, &op_0_k, &op_0_alpha, tensor_densekernel0, &op_0_ldb, tensor_dense_input, &op_0_lda, &op_0_beta, tensor_denseDense, &op_0_n);
//------ RELU
for (int id = 0; id < 64 ; id++){
tensor_denseRelu0[id] = ((tensor_denseDense[id] > 0 )? tensor_denseDense[id] : 0);
}
//--------- Gemm
char op_2_transA = 'n';
char op_2_transB = 'n';
int op_2_m = 1;
int op_2_n = 64;
int op_2_k = 64;
float op_2_alpha = 1;
float op_2_beta = 1;
int op_2_lda = 64;
int op_2_ldb = 64;
std::copy(tensor_dense_1bias0bcast, tensor_dense_1bias0bcast + 64, tensor_dense_1Dense);
BLAS::sgemm_(&op_2_transB, &op_2_transA, &op_2_n, &op_2_m, &op_2_k, &op_2_alpha, tensor_dense_1kernel0, &op_2_ldb, tensor_denseRelu0, &op_2_lda, &op_2_beta, tensor_dense_1Dense, &op_2_n);
//------ RELU
for (int id = 0; id < 64 ; id++){
tensor_dense_1Relu0[id] = ((tensor_dense_1Dense[id] > 0 )? tensor_dense_1Dense[id] : 0);
}
//--------- Gemm
char op_4_transA = 'n';
char op_4_transB = 'n';
int op_4_m = 1;
int op_4_n = 64;
int op_4_k = 64;
float op_4_alpha = 1;
float op_4_beta = 1;
int op_4_lda = 64;
int op_4_ldb = 64;
std::copy(tensor_dense_2bias0bcast, tensor_dense_2bias0bcast + 64, tensor_dense_2Dense);
BLAS::sgemm_(&op_4_transB, &op_4_transA, &op_4_n, &op_4_m, &op_4_k, &op_4_alpha, tensor_dense_2kernel0, &op_4_ldb, tensor_dense_1Relu0, &op_4_lda, &op_4_beta, tensor_dense_2Dense, &op_4_n);
//------ RELU
for (int id = 0; id < 64 ; id++){
tensor_dense_2Relu0[id] = ((tensor_dense_2Dense[id] > 0 )? tensor_dense_2Dense[id] : 0);
}
//--------- Gemm
char op_6_transA = 'n';
char op_6_transB = 'n';
int op_6_m = 1;
int op_6_n = 64;
int op_6_k = 64;
float op_6_alpha = 1;
float op_6_beta = 1;
int op_6_lda = 64;
int op_6_ldb = 64;
std::copy(tensor_dense_3bias0bcast, tensor_dense_3bias0bcast + 64, tensor_dense_3Dense);
BLAS::sgemm_(&op_6_transB, &op_6_transA, &op_6_n, &op_6_m, &op_6_k, &op_6_alpha, tensor_dense_3kernel0, &op_6_ldb, tensor_dense_2Relu0, &op_6_lda, &op_6_beta, tensor_dense_3Dense, &op_6_n);
//------ RELU
for (int id = 0; id < 64 ; id++){
tensor_dense_3Relu0[id] = ((tensor_dense_3Dense[id] > 0 )? tensor_dense_3Dense[id] : 0);
}
//--------- Gemm
char op_8_transA = 'n';
char op_8_transB = 'n';
int op_8_m = 1;
int op_8_n = 2;
int op_8_k = 64;
float op_8_alpha = 1;
float op_8_beta = 1;
int op_8_lda = 64;
int op_8_ldb = 2;
std::copy(tensor_dense_4bias0bcast, tensor_dense_4bias0bcast + 2, tensor_dense_4Dense);
BLAS::sgemm_(&op_8_transB, &op_8_transA, &op_8_n, &op_8_m, &op_8_k, &op_8_alpha, tensor_dense_4kernel0, &op_8_ldb, tensor_dense_3Relu0, &op_8_lda, &op_8_beta, tensor_dense_4Dense, &op_8_n);
for (int id = 0; id < 2 ; id++){
tensor_dense_4Sigmoid0[id] = 1 / (1 + std::exp( - tensor_dense_4Dense[id]));
}
return fTensor_dense_4Sigmoid0;
}
}; // end of Session
} //TMVA_SOFIE_Higgs_trained_model
#endif // ROOT_TMVA_SOFIE_HIGGS_TRAINED_MODEL
TF/Keras Version: 2.13.0
Generating inference code for the Keras model from Higgs_trained_model.h5 in the header Higgs_trained_model.hxx
compiling SOFIE model Higgs_trained_model
Generating inference code for the Keras model from Higgs_trained_model.h5 in the header Higgs_trained_model.hxx
size of data 10000
Number of signal entries 10000.0
Number of background entries 10000.0
Author
Lorenzo Moneta

Definition in file TMVA_SOFIE_Inference.py.