13.5.20. Macro “sensitivitySobolRankFunctionFlowrate.C”
13.5.20.1. Objective
The objective of this macro is to perform a sensitivity analysis using the SobolRank method on a set
of eight parameters used in the flowrateModel model described in
Presentation of the problem.
13.5.20.2. Macro Uranie
void sensitivitySobolRankFunctionFlowrate(const string& figure="figure.png", const string& style="",
const string& filename="", const string& treename="", int seed=0)
{
gROOT->LoadMacro("UserFunctions.C");
// Define the DataServer
TDataServer *tds = new TDataServer("tdsflowreate", "DataBase flowreate");
tds->addAttribute( new TUniformDistribution("rw", 0.05, 0.15));
tds->addAttribute( new TUniformDistribution("r", 100.0, 50000.0));
tds->addAttribute( new TUniformDistribution("tu", 63070.0, 115600.0));
tds->addAttribute( new TUniformDistribution("tl", 63.1, 116.0));
tds->addAttribute( new TUniformDistribution("hu", 990.0, 1110.0));
tds->addAttribute( new TUniformDistribution("hl", 700.0, 820.0));
tds->addAttribute( new TUniformDistribution("l", 1120.0, 1680.0));
tds->addAttribute( new TUniformDistribution("kw", 9855.0, 12045.0));
// Generation of the sample (it can be a given sample).
Int_t nS = 500;
TSampling *sampling = new TSampling(tds, "lhs", nS);
sampling->generateSample();
TLauncherFunction * tlf = new TLauncherFunction(tds, "flowrateModel");
tlf->setDrawProgressBar(kFALSE);
tlf->run();
// Create a TSobolRank object, compute indexes and print results
TSobolRank * tsobolrank = new TSobolRank(tds, "rw:r:tu:tl:hu:hl:l:kw","flowrateModel");
tsobolrank->computeIndexes("quiet");
tsobolrank->getResultTuple()->SetScanField(60);
tsobolrank->getResultTuple()->Scan("Out:Inp:Method:Order:Value:CILower:CIUpper", "","colsize=5 col=6:8::9:8:8:8");
// Print Sobol indexes
TCanvas *can = new TCanvas("c1", "Graph for the Macro sensitivitySobolRankFunctionFlowrate",5,64,1270,667);
tsobolrank->drawIndexes("Flowrate", "", "hist,first,nonewcanv");
}
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 DataServer
TDataServer *tds = new TDataServer("tdsflowreate", "DataBase flowreate");
tds->addAttribute( new TUniformDistribution("rw", 0.05, 0.15));
tds->addAttribute( new TUniformDistribution("r", 100.0, 50000.0));
tds->addAttribute( new TUniformDistribution("tu", 63070.0, 115600.0));
tds->addAttribute( new TUniformDistribution("tl", 63.1, 116.0));
tds->addAttribute( new TUniformDistribution("hu", 990.0, 1110.0));
tds->addAttribute( new TUniformDistribution("hl", 700.0, 820.0));
tds->addAttribute( new TUniformDistribution("l", 1120.0, 1680.0));
tds->addAttribute( new TUniformDistribution("kw", 9855.0, 12045.0));
SobolRank does not need a specific DOE, it works with a given sample. We generate a sample with
TSampling and we evaluate it using TLauncherFunction
Int_t nS = 500;
TSampling *sampling = new TSampling(tds, "lhs", nS);
sampling->generateSample();
TLauncherFunction * tlf = new TLauncherFunction(tds, "flowrateModel");
tlf->setDrawProgressBar(kFALSE);
tlf->run();
To instantiate the TSobolRank object, one uses the TDataServer, the name of the input and output
variables:
TSobolRank * tsobolrank = new TSobolRank(tds, "rw:r:tu:tl:hu:hl:l:kw","flowrateModel");
Computation of sensitivity indexes:
tsobolrank->computeIndexes("quiet");
It will produce one plot containing the first Sobol indices. The console is also shown below for completness.
13.5.20.3. Graph
Figure 13.36 Graph of the macro “sensitivitySobolRankFunctionFlowrate.C”
13.5.20.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
*************************************************************************************
* Row * Out * Inp * Metho * Order * Value * CILower * CIUpper *
*************************************************************************************
* 0 * flowra * rw * Sobol * First * 0.816872 * -1 * -1 *
* 1 * flowra * r * Sobol * First * -0.04759 * -1 * -1 *
* 2 * flowra * tu * Sobol * First * -0.04164 * -1 * -1 *
* 3 * flowra * tl * Sobol * First * -0.00414 * -1 * -1 *
* 4 * flowra * hu * Sobol * First * 0.041272 * -1 * -1 *
* 5 * flowra * hl * Sobol * First * 0.036748 * -1 * -1 *
* 6 * flowra * l * Sobol * First * 0.062867 * -1 * -1 *
* 7 * flowra * kw * Sobol * First * 0.021273 * -1 * -1 *
*************************************************************************************