13.8.1. Macro “reoptimizeHollowBarCode.C

13.8.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.8.1.2. Macro Uranie

    using namespace URANIE::DataServer;
    using namespace URANIE::Relauncher;
    using namespace URANIE::Reoptimizer;

    // variables
    TAttribute x("x", 0.0, 1.0),
      y("y", 0.0, 1.0),
      thick("thick"), // thickness
      sect("sect"), // section of the pipe
      dist("dist"); // distortion

    // Creating the TCodeEval, dumping output of the dummy python in an output file
    string python_exec = "python3";
    if(string(gSystem->GetBuildArch()) == "win64")
        python_exec.pop_back();
    TCodeEval code((python_exec +" bar.py > bartoto.dat").data());

    // Pass the python script itself as a input file. x and y will be modified in bar.py directly
    TKeyScript inputfile("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 variables
    TFlatResult outputfile("bartoto.dat");
    outputfile.addOutput(&thick);
    outputfile.addOutput(&sect);
    outputfile.addOutput(&dist);
    code.addOutputFile(&outputfile);

    // Create a runner
    TSequentialRun runner(&code);
    runner.startSlave(); // Usual Relauncher construction

    if(runner.onMaster())
    {
        // Create the TDS
        TDataServer tds("vizirDemo", "Param de l'opt vizir pour la barre");
        tds.addAttribute(&x);
        tds.addAttribute(&y);

        // Choose a solver
        TNloptCobyla solv;
        //TNloptBobyqa solv;
        //TNloptPraxis solv;
        //TNloptNelderMead solv;
        ///TNloptSubplexe solv;

        // Create the single-objective constrained optimizer
        TNlopt opt(&tds, &runner, &solv);

        // add the objective
        opt.addObjective(&sect); // minimizing the section

	    // and the constrains
        TLesserFit constrDist(14);
        opt.addConstraint(&dist,&constrDist); // on the distortion (dist < 14)
        TGreaterFit positiv(0.4);
        opt.addConstraint(&thick,&positiv); // and on the thickness (thick > 0.4)

        // Starting point and maximum evaluation
        vector<double> point{0.9 , 0.2};
        opt.setStartingPoint(point.size(),&point[0]);
        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
    TAttribute x("x", 0.0, 1.0),
      y("y", 0.0, 1.0),
      thick("thick"), // thickness
      sect("sect"), // section of the pipe
      dist("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 an output file
    TCodeEval code((python_exec +" bar.py > bartoto.dat").data());

    // Pass the python script itself as a input file. x and y will be modified in bar.py directly
    TKeyScript inputfile("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 variables
    TFlatResult outputfile("bartoto.dat");
    outputfile.addOutput(&thick);
    outputfile.addOutput(&sect);
    outputfile.addOutput(&dist);
    code.addOutputFile(&outputfile);

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

        // Choose a solver
        TNloptCobyla solv;
        //TNloptBobyqa solv;
        //TNloptPraxis solv;
        //TNloptNelderMead solv;
        ///TNloptSubplexe solv;

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
        TNlopt opt(&tds, &runner, &solv);

        // add the objective
        opt.addObjective(&sect); // minimizing the section

	    // and the constrains
        TLesserFit constrDist(14);
        opt.addConstraint(&dist,&constrDist); // on the distortion (dist < 14)
        TGreaterFit positiv(0.4);
        opt.addConstraint(&thick,&positiv); // and on the 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
        vector<double> point{0.9 , 0.2};
        opt.setStartingPoint(point.size(),&point[0]);
        opt.setMaximumEval(1000);

        opt.solverLoop(); // running the optimization

13.8.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 *
***************************************************************************