13.7.5. Macro “modelerFlowrateNeuralNetworksLoadingPMML.py

13.7.5.1. Objective

The objective of this macro is to build a surrogate model (an Artificial Neural Network) from a PMML file.

13.7.5.2. Macro Uranie

"""
Example of neural network loaded from pmml for flowrate
"""
from URANIE import DataServer, Launcher, Modeler
import ROOT

# Create a DataServer.TDataServer
# Load a database in an ASCII file
tds = DataServer.TDataServer("tdsFlowrate", "tds for flowrate")
tds.fileDataRead("_flowrate_sampler_launcher_app_.dat")

c = ROOT.TCanvas("c1", "Graph for the Macro modeler", 5, 64, 1270, 667)
# Build a surrogate model (Artificial Neural Networks) from the PMML file
tann = Modeler.TANNModeler(tds, "uranie_ann_flowrate.pmml", "ANNflowrate")

# export the surrogate model in a C file
tann.exportFunction("c++", "uranie_ann_flowrate_loaded", "ANNflowrate")
# load the surrogate model in the C file
ROOT.gROOT.LoadMacro("uranie_ann_flowrate_loaded.C")

tdsv = DataServer.TDataServer("tdsv", "tds for surrogate model")
tdsv.fileDataRead("_flowrate_sampler_launcher_val_.dat")

print(tdsv.getNPatterns())

# evaluate the surrogate model on the database
tlf = Launcher.TLauncherFunction(tdsv, "ANNflowrate",
                                 "rw:r:tu:tl:hu:hl:l:kw", "yann")
tlf.run()

tdsv.startViewer()

tdsv.Draw("yann:yhat")
# tdsv.draw("yhat")

A TDataServer loads _flowrate_sampler_launcher_app_.dat:

tds = DataServer.TDataServer("tdsFlowrate", "tds for flowrate")
tds.fileDataRead("_flowrate_sampler_launcher_app_.dat")

The surrogate model is loaded from a PMML file:

tann = Modeler.TANNModeler(tds, "uranie_ann_flowrate.pmml", "ANNflowrate")

The model is exported in an external file in C++ language "uranie_ann_flowrate_loaded.C" where the function name is ANNflowrate:

tann.exportFunction("c++", "uranie_ann_flowrate_loaded", "ANNflowrate")

The model is loaded from the macro "uranie_ann_flowrate_loaded.C" and applied on the database with the function ANNflowrate:

ROOT.gROOT.LoadMacro("uranie_ann_flowrate_loaded.C")
tdsv = DataServer.TDataServer("tdsv", "tds for surrogate model")
tdsv.fileDataRead("_flowrate_sampler_launcher_val_.dat")
tlf = Launcher.TLauncherFunction(tdsv, "ANNflowrate",
                                 "rw:r:tu:tl:hu:hl:l:kw", "yann")
tlf.run()

13.7.5.3. Graph

../../_images/modelerFlowrateNeuralNetworksLoadingPMML.png

Figure 13.41 Graph of the macro “modelerFlowrateNeuralNetworksLoadingPMML.py”

13.7.5.4. Console

--- Uranie v4.11/0 --- Developed with ROOT (6.36.06)
                      Copyright (C) 2013-2026 CEA/DES 
                      Contact: support-uranie@cea.fr 
                      Date: Thu Feb 12, 2026

PMML Constructor: uranie_ann_flowrate.pmml

*******************************
*** TModeler::exportFunction lang[c++] file[uranie_ann_flowrate_loaded] name[ANNflowrate] soption[]

*******************************

*******************************
*** exportFunction lang[c++] file[uranie_ann_flowrate_loaded] name[ANNflowrate]
*** End Of exportFunction
*******************************