13.5.16. Macro “sensitivityJohnsonRWFunctionFlowrate.C

13.5.16.1. Objective

The objective of this macro is to perform a sensitivity analysis using the Johnson’s relative weight method on a set of eight parameters used in the flowrateModel model described in Presentation of the problem.

13.5.16.2. Macro Uranie

    // Define the DataServer
    TDataServer *tds = new TDataServer("tdsflowrate", "DataBase flowrate");

    gROOT->LoadMacro("UserFunctions.C");
          
    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 = 1000;
    TString FuncName="flowrateModel";

    TJohnsonRW *tjrw = new TJohnsonRW(tds, FuncName, nS, "rw:r:tu:tl:hu:hl:l:kw", FuncName);
    tjrw->computeIndexes();

    // Get the results on screen 
    tjrw->getResultTuple()->Scan("Out:Inp:Method:Value","Order==\"First\"");

    // Get the results as plots
    TCanvas *cc = new TCanvas("canhist", "histgramme");
    tjrw->drawIndexes("Flowrate", "", "nonewcanv,hist,first"); 
    cc->Print("appliUranieFlowrateJohnsonRW1000Histogram.png"); 

    TCanvas *ccc = new TCanvas("canpie", "TPie");
    tjrw->drawIndexes("Flowrate", "", "nonewcanv,pie"); 
    ccc->Print("appliUranieFlowrateJohnsonRW1000Pie.png");

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("tdsflowrate", "DataBase flowrate");
    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 TJohnsonRW object, one uses the TDataServer, the name of the function, the number of samplings needed to perform sensitivity analysis (here \(n_S\) =4000) and the input and output variables:

    TJohnsonRW *tjrw = new TJohnsonRW(tds, FuncName, nS, "rw:r:tu:tl:hu:hl:l:kw", FuncName);

Computation of sensitivity indexes:

    tjrw->computeIndexes();

The rest is very common to all sensitivity macros discussed here: it will produce two plots (the first one being a histogram show below) and the console is also shown below for completness.

13.5.16.3. Graph

../../_images/appliUranieFlowrateJohnsonRW1000Histogram.png

Figure 13.32 Graph of the macro “sensitivityJohnsonRWFunctionFlowrate.C”

13.5.16.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 [tdsflowrate] contains data: we need to empty it ! 
 <URANIE::INFO> *** END of URANIE INFORMATION ***
 <URANIE::INFO> 
************************************************************
*    Row   *       Out *       Inp *    Method *     Value *
************************************************************
*        0 * flowrateM *        rw * JohnsonRW * 0.8636266 *
*        2 * flowrateM *         r * JohnsonRW * 4.598e-05 *
*        4 * flowrateM *        tu * JohnsonRW * 8.775e-06 *
*        6 * flowrateM *        tl * JohnsonRW * 9.689e-05 *
*        8 * flowrateM *        hu * JohnsonRW * 0.0439173 *
*       10 * flowrateM *        hl * JohnsonRW * 0.0435594 *
*       12 * flowrateM *         l * JohnsonRW * 0.0378304 *
*       14 * flowrateM *        kw * JohnsonRW * 0.0109144 *
*       16 * flowrateM *    __R2__ * JohnsonRW * 0.9447092 *
*       18 * flowrateM *   __R2A__ * JohnsonRW * 0.9442633 *
************************************************************
==> 10 selected entries