13.5.4. Macro “sensitivityFASTFunctionFlowrate.C”
13.5.4.1. Objective
The objective of this macro is to perform a Fast sensitivity analysis on a set of eight parameters
used in the flowrateModel model described in Presentation of the problem.
13.5.4.2. Macro Uranie
void sensitivityFASTFunctionFlowrate(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));
// \param Size of a sampling.
Int_t nS = 4000;
// Graph
TFast * tfast = new TFast(tds, "flowrateModel", nS);
tfast->setDrawProgressBar(kFALSE);
tfast->computeIndexes("graph");
tfast->getResultTuple()->Scan("Out:Inp:Order:Method:Value","Algo==\"--first--\"");
}
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:
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));
To instantiate the TFast object, one uses the TDataServer, the name of the function and the number of
samplings needed to perform sensitivity analysis (here nS=500):
TFast * tfast = new TFast(tds, "flowrateModel", nS);
Computation of sensitivity indexes:
tfast->computeIndexes("graph");
13.5.4.3. Graph
Figure 13.20 Graph of the macro “sensitivityFASTFunctionFlowrate.C”
13.5.4.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::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/TSpaceFilling.cxx] Line[167]
<URANIE::INFO> TSamplerStochastic::init: the TDS [tdsflowreate] contains data: we need to empty it !
<URANIE::INFO> *** END of URANIE INFORMATION ***
<URANIE::INFO>
************************************************************************
* Row * Out * Inp * Order * Method * Value *
************************************************************************
* 0 * flowrateM * rw * First * FAST * 0.8278187 *
* 2 * flowrateM * r * First * FAST * 8.924e-07 *
* 4 * flowrateM * tu * First * FAST * 2.308e-06 *
* 6 * flowrateM * tl * First * FAST * 3.204e-05 *
* 8 * flowrateM * hu * First * FAST * 0.0414390 *
* 10 * flowrateM * hl * First * FAST * 0.0414046 *
* 12 * flowrateM * l * First * FAST * 0.0392873 *
* 14 * flowrateM * kw * First * FAST * 0.0094983 *
************************************************************************
==> 8 selected entries