13.5.10. Macro “sensitivitySobolFunctionFlowrateRunner.C”
13.5.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.5.10.2. Macro Uranie
gROOT->LoadMacro("UserFunctions.C");
// Define the attributes
TUniformDistribution rw("rw", 0.05, 0.15);
TUniformDistribution r("r", 100.0, 50000.0);
TUniformDistribution tu("tu", 63070.0, 115600.0);
TUniformDistribution tl("tl", 63.1, 116.0);
TUniformDistribution hu("hu", 990.0, 1110.0);
TUniformDistribution hl("hl", 700.0, 820.0);
TUniformDistribution l("l", 1120.0, 1680.0);
TUniformDistribution kw("kw", 9855.0, 12045.0);
// Create the evaluator
TCIntEval code("flowrateModel");
// Create output attribute
TAttribute yout("flowrateModel");
// Provide input/output attributes to the assessor
code.setInputs(8, &rw, &r, &tu, &tl, &hu, &hl, &l, &kw);
code.setOutputs(1, &yout);
TSequentialRun run(&code); // To be replaced to distribute the computation
run.startSlave();
if( run.onMaster())
{
// Create the dataserver
TDataServer *tds = new TDataServer("sobol", "foo bar pouet chocolat");
tds->addAttribute(&rw);
tds->addAttribute(&r);
tds->addAttribute(&tu);
tds->addAttribute(&tl);
tds->addAttribute(&hu);
tds->addAttribute(&hl);
tds->addAttribute(&l);
tds->addAttribute(&kw);
// Create the sobol object
Int_t ns = 100000;
TSobol * tsobol = new TSobol(tds, &run, ns);
tsobol->setDrawProgressBar(kFALSE);
tsobol->computeIndexes();
TCanvas *cc = new TCanvas("c1", "histgramme",5,64,1270,667);
TPad *pad = new 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)
gROOT->LoadMacro("UserFunctions.C");
Each parameter is related to the TDataServer as a TAttribute and obeys an uniform law on specific interval:
// Define the attributes
TUniformDistribution rw("rw", 0.05, 0.15);
TUniformDistribution r("r", 100.0, 50000.0);
TUniformDistribution tu("tu", 63070.0, 115600.0);
TUniformDistribution tl("tl", 63.1, 116.0);
TUniformDistribution hu("hu", 990.0, 1110.0);
TUniformDistribution hl("hl", 700.0, 820.0);
TUniformDistribution l("l", 1120.0, 1680.0);
TUniformDistribution kw("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
TCIntEval code("flowrateModel");
// Create output attribute
TAttribute yout("flowrateModel");
// Provide input/output attributes to the assessor
code.setInputs(8, &rw, &r, &tu, &tl, &hu, &hl, &l, &kw);
code.setOutputs(1, &yout);
TSequentialRun run(&code); // To be replaced to distribute the 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 * tsobol = new 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.5.10.3. Graph
Figure 13.26 Graph of the macro “sensitivitySobolFunctionFlowrateRunner.C”
13.5.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]