13.6.10. Macro “sensitivitySobolFunctionFlowrateRunner.py

13.6.10.1. Objective

The objective of this macro is to perform Sobol sensitivity analysis on a set of eight parameters used in the flowrateModel model described in Presentation of the problem, but this time using the Relauncher architecture.

13.6.10.2. Macro Uranie

"""
Example of Sobol estimation for the flowrate function with Relauncher approach
"""
from URANIE import DataServer, Relauncher, Sensitivity
import ROOT

ROOT.gROOT.LoadMacro("UserFunctions.C")

# Define the DataServer
rw = DataServer.TUniformDistribution("rw", 0.05, 0.15)
r = DataServer.TUniformDistribution("r", 100.0, 50000.0)
tu = DataServer.TUniformDistribution("tu", 63070.0, 115600.0)
tl = DataServer.TUniformDistribution("tl", 63.1, 116.0)
hu = DataServer.TUniformDistribution("hu", 990.0, 1110.0)
hl = DataServer.TUniformDistribution("hl", 700.0, 820.0)
lvar = DataServer.TUniformDistribution("l", 1120.0, 1680.0)
kw = DataServer.TUniformDistribution("kw", 9855.0, 12045.0)

# Create the evaluator
code = Relauncher.TCIntEval("flowrateModel")
# Create output attribute
yout = DataServer.TAttribute("flowrateModel")
# Provide input/output attributes to the assessor
code.addInput(rw)
code.addInput(r)
code.addInput(tu)
code.addInput(tl)
code.addInput(hu)
code.addInput(hl)
code.addInput(lvar)
code.addInput(kw)
code.addOutput(yout)

run = Relauncher.TSequentialRun(code)  # Replace to distribute computation
run.startSlave()
if run.onMaster():
    # Create the dataserver
    tds = DataServer.TDataServer("sobol", "foo bar pouet chocolat")
    code.addAllInputs(tds)

    # Create the sobol object
    ns = 100000
    tsobol = Sensitivity.TSobol(tds, run, ns)
    tsobol.setDrawProgressBar(ROOT.kFALSE)
    tsobol.computeIndexes()

    cc = ROOT.TCanvas("c1", "histgramme", 5, 64, 1270, 667)
    pad = ROOT.TPad("pad", "pad", 0, 0.03, 1, 1)
    pad.Draw()
    pad.Divide(2, 1)
    pad.cd(1)
    tsobol.drawIndexes("Flowrate", "", "nonewcanv, hist, all")

    pad.cd(2)
    tsobol.drawIndexes("Flowrate", "", "nonewcanv, pie, first")

    tds.exportData("_sobol_launching_.dat")

The function flowrateModel is loaded from the macro UserFunctions.C (the file can be found in ${URANIESYS}/share/uranie/macros)

ROOT.gROOT.LoadMacro("UserFunctions.C")

Each parameter is related to the TDataServer as a TAttribute and obeys an uniform law on specific interval:

# Define the DataServer
rw = DataServer.TUniformDistribution("rw", 0.05, 0.15)
r = DataServer.TUniformDistribution("r", 100.0, 50000.0)
tu = DataServer.TUniformDistribution("tu", 63070.0, 115600.0)
tl = DataServer.TUniformDistribution("tl", 63.1, 116.0)
hu = DataServer.TUniformDistribution("hu", 990.0, 1110.0)
hl = DataServer.TUniformDistribution("hl", 700.0, 820.0)
lvar = DataServer.TUniformDistribution("l", 1120.0, 1680.0)
kw = DataServer.TUniformDistribution("kw", 9855.0, 12045.0)

The interface to the function is then defined, using the Relauncher interface, through a TCIntEval object and a sequential runner:

# Create the evaluator
code = Relauncher.TCIntEval("flowrateModel")
# Create output attribute
yout = DataServer.TAttribute("flowrateModel")
# Provide input/output attributes to the assessor
code.addInput(rw)
code.addInput(r)
code.addInput(tu)
code.addInput(tl)
code.addInput(hu)
code.addInput(hl)
code.addInput(lvar)
code.addInput(kw)
code.addOutput(yout)

run = Relauncher.TSequentialRun(code)  # Replace to distribute computation
run.startSlave()

To instantiate the TSobol, one uses the TDataServer, a pointer to the runner and the number of samplings needed to perform sensitivity analysis (here ns=600):

tsobol = Sensitivity.TSobol(tds, run, ns)

Computation of the sensitivity indexes:

tsobol.computeIndexes()

Data are exported from the TDataServer to an ASCII file:

tds.exportData("_sobol_launching_.dat")

13.6.10.3. Graph

../../_images/sensitivitySobolFunctionFlowrateRunner.png

Figure 13.26 Graph of the macro “sensitivitySobolFunctionFlowrateRunner.py”

13.6.10.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

 <URANIE::WARNING> 
 <URANIE::WARNING> *** URANIE WARNING ***
 <URANIE::WARNING> *** File[${SOURCEDIR}/dataSERVER/souRCE/TDataServer.cxx] Line[8531]
 <URANIE::WARNING> TDataServer::getTuple Error : There is no tree!
 <URANIE::WARNING> *** END of URANIE WARNING ***
 <URANIE::WARNING> 
 <URANIE::INFO> 
 <URANIE::INFO> *** URANIE INFORMATION ***
 <URANIE::INFO> *** File[${SOURCEDIR}/meTIER/sampler/souRCE/TSamplerStochastic.cxx] Line[66]
 <URANIE::INFO> TSamplerStochastic::init: the TDS [sobol] contains data: we need to empty it ! 
 <URANIE::INFO> *** END of URANIE INFORMATION ***
 <URANIE::INFO> 
 ** Case of Output atty [flowrateModel]  nSimPerIndex 10000
 ** Input att [rw] First [0.830033] Total Order[0.865762]
 ** Input att [r] First [0] Total Order[0.000102212]
 ** Input att [tu] First [0] Total Order[0.0001]
 ** Input att [tl] First [0] Total Order[0.000110756]
 ** Input att [hu] First [0.0417298] Total Order[0.0554922]
 ** Input att [hl] First [0.0367345] Total Order[0.0526188]
 ** Input att [l] First [0.0384214] Total Order[0.0535728]
 ** Input att [kw] First [0.00669831] Total Order[0.0132316]