13.5.19. Macro “sensitivityHSICFunctionFlowrate.C

13.5.19.1. Objective

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

13.5.19.2. Macro Uranie

void sensitivityHSICFunctionFlowrate(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 THSIC object, compute indexes and print results
    THSIC * thsic = new THSIC(tds, "rw:r:tu:tl:hu:hl:l:kw","flowrateModel");
    thsic->computeIndexes("quiet");
    thsic->getResultTuple()->SetScanField(60);
    thsic->getResultTuple()->Scan("Out:Inp:Method:Order:Value:CILower:CIUpper", "","colsize=5 col=6:8::9:8:8:8");

    // Print HSIC indexes
    TCanvas  *can = new TCanvas("c1", "Graph for the Macro sensitivityHSICFunctionFlowrate",5,64,1270,667);
    thsic->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));

HSIC 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 THSIC object, one uses the TDataServer, the name of the input and output variables:

    THSIC * thsic = new THSIC(tds, "rw:r:tu:tl:hu:hl:l:kw","flowrateModel");

Computation of sensitivity indexes:

    thsic->computeIndexes("quiet");

It will produce one plots containing the HSIC indexes and the p-value to test the Independance between inputs and outputs.The console is also shown below for completness.

13.5.19.3. Graph

../../_images/sensitivityHSICFunctionFlowrate.png

Figure 13.35 Graph of the macro “sensitivityHSICFunctionFlowrate.C”

13.5.19.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 *  HSic *    R2HSic * 0.804431 *       -1 *       -1 *
*        1 * flowra *       rw *  HSic *      HSic * 0.072723 *       -1 *       -1 *
*        2 * flowra *       rw *  HSic *   pValues * 6.8e-106 *       -1 *       -1 *
*        3 * flowra *        r *  HSic *    R2HSic * 0.002238 *       -1 *       -1 *
*        4 * flowra *        r *  HSic *      HSic * 0.000202 *       -1 *       -1 *
*        5 * flowra *        r *  HSic *   pValues * 0.712174 *       -1 *       -1 *
*        6 * flowra *       tu *  HSic *    R2HSic * 0.002528 *       -1 *       -1 *
*        7 * flowra *       tu *  HSic *      HSic * 0.000228 *       -1 *       -1 *
*        8 * flowra *       tu *  HSic *   pValues * 0.662668 *       -1 *       -1 *
*        9 * flowra *       tl *  HSic *    R2HSic * 0.003806 *       -1 *       -1 *
*       10 * flowra *       tl *  HSic *      HSic * 0.000344 *       -1 *       -1 *
*       11 * flowra *       tl *  HSic *   pValues * 0.449068 *       -1 *       -1 *
*       12 * flowra *       hu *  HSic *    R2HSic * 0.025554 *       -1 *       -1 *
*       13 * flowra *       hu *  HSic *      HSic * 0.002309 *       -1 *       -1 *
*       14 * flowra *       hu *  HSic *   pValues * 3.13e-05 *       -1 *       -1 *
*       15 * flowra *       hl *  HSic *    R2HSic * 0.020243 *       -1 *       -1 *
*       16 * flowra *       hl *  HSic *      HSic * 0.001829 *       -1 *       -1 *
*       17 * flowra *       hl *  HSic *   pValues * 0.000445 *       -1 *       -1 *
*       18 * flowra *        l *  HSic *    R2HSic * 0.024934 *       -1 *       -1 *
*       19 * flowra *        l *  HSic *      HSic * 0.002253 *       -1 *       -1 *
*       20 * flowra *        l *  HSic *   pValues * 5.64e-05 *       -1 *       -1 *
*       21 * flowra *       kw *  HSic *    R2HSic * 0.010567 *       -1 *       -1 *
*       22 * flowra *       kw *  HSic *      HSic * 0.000955 *       -1 *       -1 *
*       23 * flowra *       kw *  HSic *   pValues * 0.032459 *       -1 *       -1 *
*************************************************************************************