13.9.1. Macro “reoptimizeHollowBarCode.py

13.9.1.1. Objective

The objective of the macro is to optimize the section of the hollow bar defined in Sizing of a hollow bar example problem using the NLopt solvers (reducing it to a single-criterion optimisation as already explained in Local solver. This can be done with different solvers, the results being achieved within more or less time and following the requested constraints with more or less accuracy (depending on the hypothesis embedded by the chosen solver).

13.9.1.2. Macro Uranie

"""
Example of hollow bar optimisation
"""
import numpy as np
from URANIE import DataServer, Relauncher, Reoptimizer
import ROOT

# variables
x = DataServer.TAttribute("x", 0.0, 1.0)
y = DataServer.TAttribute("y", 0.0, 1.0)
thick = DataServer.TAttribute("thick")  # thickness
sect = DataServer.TAttribute("sect")  # section of the pipe
dist = DataServer.TAttribute("dist")  # distortion

# Creating the TCodeEval, dumping output of the dummy python in output file
python_exec = "python3"
if ROOT.gSystem.GetBuildArch() == "win64":
    python_exec = python_exec[:-1]
code = Relauncher.TCodeEval(python_exec + " bar.py > bartoto.dat")

# Pass the python script itself as input. inputs are modified in bar.py
inputfile = Relauncher.TKeyScript("bar.py")
inputfile.addInput(x, "x")
inputfile.addInput(y, "y")
code.addInputFile(inputfile)

# precise the name of the output file in which to read the three output
outputfile = Relauncher.TFlatResult("bartoto.dat")
outputfile.addOutput(thick)
outputfile.addOutput(sect)
outputfile.addOutput(dist)
code.addOutputFile(outputfile)

# Create a runner
runner = Relauncher.TSequentialRun(code)
runner.startSlave()  # Usual Relauncher construction

if runner.onMaster():

    # Create the TDS
    tds = DataServer.TDataServer("vizirDemo", "Param de l'opt vizir")
    tds.addAttribute(x)
    tds.addAttribute(y)

    # Choose a solver
    solv = Reoptimizer.TNloptCobyla()
    # solv = Reoptimizer.TNloptBobyqa()
    # solv = Reoptimizer.TNloptPraxis()
    # solv = Reoptimizer.TNloptNelderMead()
    # solv = Reoptimizer.TNloptSubplexe()

    # Create the single-objective constrained optimizer
    opt = Reoptimizer.TNlopt(tds, runner, solv)

    # add the objective
    opt.addObjective(sect)  # minimizing the section

    # and the constrains
    constrDist = Reoptimizer.TLesserFit(14)
    opt.addConstraint(dist, constrDist)  # distortion (dist < 14)
    positiv = Reoptimizer.TGreaterFit(0.4)
    opt.addConstraint(thick, positiv)  # thickness (thick > 0.4)

    # Starting point and maximum evaluation
    point = np.array([0.9, 0.2])
    opt.setStartingPoint(len(point), point)
    opt.setMaximumEval(1000)

    opt.solverLoop()  # running the optimization

    # Stop the slave processes
    runner.stopSlave()

# solution
tds.getTuple().Scan("*", "", "colsize=9 col=:::5:4")

The variables are defined as follow:

# variables
x = DataServer.TAttribute("x", 0.0, 1.0)
y = DataServer.TAttribute("y", 0.0, 1.0)
thick = DataServer.TAttribute("thick")  # thickness
sect = DataServer.TAttribute("sect")  # section of the pipe
dist = DataServer.TAttribute("dist")  # distortion

where the first two are the input ones while the last ones are computed using the provided code (as explained in Sizing of a hollow bar example problem). This code is configured through these lines:

# Creating the TCodeEval, dumping output of the dummy python in output file
code = Relauncher.TCodeEval(python_exec + " bar.py > bartoto.dat")

# Pass the python script itself as input. inputs are modified in bar.py
inputfile = Relauncher.TKeyScript("bar.py")
inputfile.addInput(x, "x")
inputfile.addInput(y, "y")
code.addInputFile(inputfile)

# precise the name of the output file in which to read the three output
outputfile = Relauncher.TFlatResult("bartoto.dat")
outputfile.addOutput(thick)
outputfile.addOutput(sect)
outputfile.addOutput(dist)
code.addOutputFile(outputfile)

# Create a runner
runner = Relauncher.TSequentialRun(code)
runner.startSlave()  # Usual Relauncher construction

The usual Relauncher construction is followed, using a TSequentialRun runner and the solver is chosen in these lines:

# Choose a solver
solv = Reoptimizer.TNloptCobyla()
# solv = Reoptimizer.TNloptBobyqa()
# solv = Reoptimizer.TNloptPraxis()
# solv = Reoptimizer.TNloptNelderMead()
# solv = Reoptimizer.TNloptSubplexe()

Combining the runner, solver and dataserver, the master object is created and the objective and constraint are defined (keeping in mind that only single-criterion problems are implemented when dealing with NLopt, so the distortion criteria is downgraded to a constraint). This is done in

# Create the single-objective constrained optimizer
opt = Reoptimizer.TNlopt(tds, runner, solv)

# add the objective
opt.addObjective(sect)  # minimizing the section

# and the constrains
constrDist = Reoptimizer.TLesserFit(14)
opt.addConstraint(dist, constrDist)  # distortion (dist < 14)
positiv = Reoptimizer.TGreaterFit(0.4)
opt.addConstraint(thick, positiv)  # thickness (thick > 0.4)

Finally the starting point is set along with the maximal number of evaluation just before starting the loop.

# Starting point and maximum evaluation
point = np.array([0.9, 0.2])
opt.setStartingPoint(len(point), point)
opt.setMaximumEval(1000)

opt.solverLoop()  # running the optimization

13.9.1.3. Console

This macro leads to the following result

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

|....:....|....:....|....:....|....:....
|***************************************************************************
*    Row   * vizirDemo *       x.x *       y.y * thick * sect * dist.dist *
***************************************************************************
*        0 *         0 * 0.5173156 * 0.1173173 * 0.399 * 0.25 * 13.999986 *
***************************************************************************