{ "cells": [ { "cell_type": "markdown", "id": "cf3f752e", "metadata": {}, "source": [ "# TMVA_SOFIE_RSofieReader\n", "This macro provides an example of using a trained model with Keras\n", "and make inference using SOFIE with the RSofieReader class\n", "This macro uses as input a Keras model generated with the\n", "TMVA_Higgs_Classification.C tutorial\n", "You need to run that macro before to generate the trained Keras model\n", "\n", "\n", "Execute in this order:\n", "```\n", "root TMVA_Higgs_Classification.C\n", "root TMVA_SOFIE_RSofieReader.C\n", "```\n", "\n", "\n", "\n", "**Author:** Lorenzo Moneta \n", "This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, May 19, 2026 at 08:23 PM." ] }, { "cell_type": "code", "execution_count": null, "id": "4f33b656", "metadata": { "collapsed": false }, "outputs": [], "source": [ "using namespace TMVA::Experimental;" ] }, { "cell_type": "code", "execution_count": null, "id": "3d6181e7", "metadata": { "collapsed": false }, "outputs": [], "source": [ "RSofieReader model(\"HiggsModel.keras\", {}, true );" ] }, { "cell_type": "markdown", "id": "9193917d", "metadata": {}, "source": [ "for debugging\n", "RSofieReader model(\"Higgs_trained_model.keras\", {}, true);" ] }, { "cell_type": "markdown", "id": "a77f7cc0", "metadata": {}, "source": [ "the input shape for this model is a tensor with shape (1,7)" ] }, { "cell_type": "code", "execution_count": null, "id": "5e250257", "metadata": { "collapsed": false }, "outputs": [], "source": [ "std::vector input = {0.1,0.2,0.3,0.4,0.5,0.6,0.7};" ] }, { "cell_type": "markdown", "id": "3d70f2ca", "metadata": {}, "source": [ "predict model on a single event (takes a std::vector)" ] }, { "cell_type": "code", "execution_count": null, "id": "482f6871", "metadata": { "collapsed": false }, "outputs": [], "source": [ "auto output = model.Compute(input);\n", "\n", "std::cout << \"Event prediction = \" << output[0] << std::endl;" ] }, { "cell_type": "markdown", "id": "b7dc8059", "metadata": {}, "source": [ "predict model now on a input file using RDataFrame" ] }, { "cell_type": "code", "execution_count": null, "id": "5df4fdf7", "metadata": { "collapsed": false }, "outputs": [], "source": [ "std::string inputFileName = \"Higgs_data.root\";\n", "std::string inputFile = std::string{gROOT->GetTutorialDir()} + \"/machine_learning/data/\" + inputFileName;\n", "\n", "\n", "ROOT::RDataFrame df1(\"sig_tree\", inputFile);\n", "\n", "auto h1 = df1.Define(\"DNN_Values\", Compute<7, float>(model),\n", " {\"m_jj\", \"m_jjj\", \"m_lv\", \"m_jlv\", \"m_bb\", \"m_wbb\", \"m_wwbb\"})\n", " .Define(\"y\",\"DNN_Values[0]\")\n", " .Histo1D({\"h_sig\", \"\", 100, 0, 1}, \"y\");\n", "\n", "ROOT::RDataFrame df2(\"bkg_tree\", inputFile);\n", "auto h2 = df2.Define(\"DNN_Values\", Compute<7, float>(model),\n", " {\"m_jj\", \"m_jjj\", \"m_lv\", \"m_jlv\", \"m_bb\", \"m_wbb\", \"m_wwbb\"})\n", " .Define(\"y\",\"DNN_Values[0]\")\n", " .Histo1D({\"h_bkg\", \"\", 100, 0, 1}, \"y\");\n", "\n", "h1->SetLineColor(kRed);\n", "h2->SetLineColor(kBlue);\n", "\n", "auto c1 = new TCanvas();\n", "gStyle->SetOptStat(0);\n", "\n", "h2->DrawClone();\n", "h1->DrawClone(\"SAME\");\n", "c1->BuildLegend();" ] } ], "metadata": { "kernelspec": { "display_name": "ROOT C++", "language": "c++", "name": "root" }, "language_info": { "codemirror_mode": "text/x-c++src", "file_extension": ".C", "mimetype": " text/x-c++src", "name": "c++" } }, "nbformat": 4, "nbformat_minor": 5 }