13.5.7. Macro “sensitivityMorrisFunctionFlowrateRunner.C

13.5.7.1. Objective

The objective of this macro is to perform a Morris 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.7.2. Macro Uranie

void sensitivityMorrisFunctionFlowrateRunner(Int_t nk = 5, const string& figure="figure.png", 
    const string& style="", const string& filename="", const string& treename="", int seed=0)
{
    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);    
        Int_t nreplique = 3;
        Int_t nlevel = 10;
        // Create the Morris object     
        TMorris * scmo = new TMorris(tds, &run, nreplique, nlevel); 
        scmo->setDrawProgressBar(kFALSE);
        scmo->generateSample();
        
        tds->exportData("_morris_sampling_.dat");
        scmo->computeIndexes();
        
        tds->exportData("_morris_launching_.dat");

        TTree *ntresu = scmo->getMorrisResults();
        ntresu->Scan("*");
         
        // Graph
        TCanvas  *cc = new TCanvas("c1", "Graph for the Macro sensitivityMorrisFunctionFlowrateRunner",5,64,1270,667);
        TPad *pad = new TPad("pad","pad",0, 0.03, 1, 1); pad->Draw();
        pad->Divide(2);
        pad->cd(1);
        scmo->drawSample("", -1,"nonewcanv");
        pad->cd(2);
        scmo->drawIndexes("mustar,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 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();

The dataserver object is defined only on the master to avoid useless replication if one wants to run the estimation of the function in parallel (by changing the TSequentialRun by either a TThreadedRun or a TMpiRun). To instantiate the TMorris object, one uses the TDataServer, a pointer to the chosen runner, the number of replicas (here nreplique=3), the level parameter (here nlevel=10)

    TMorris * scmo = new TMorris(tds, &run, nreplique, nlevel); 

Creation of the sampling:

    scmo->generateSample();

Data are exported in an ASCII file:

    tds->exportData("_morris_sampling_.dat");

Computation of sensitivity indexes:

    scmo->computeIndexes();

The rest of the code is providing command to get a final plot.

13.5.7.3. Graph

../../_images/sensitivityMorrisFunctionFlowrateRunner.png

Figure 13.23 Graph of the macro “sensitivityMorrisFunctionFlowrateRunner.C”

13.5.7.4. Console

************************************************************************
*    Row   *     Input *    Output *     mu.mu * mustar.mu * sigma.sig *
************************************************************************
*        0 *        rw * flowrateM * 127.47900 * 127.47900 * 34.521839 *
*        1 *         r * flowrateM * -0.069601 * 0.0696013 * 0.0793689 *
*        2 *        tu * flowrateM * 0.0004201 * 0.0004201 * 0.0004641 *
*        3 *        tl * flowrateM * 0.4659763 * 0.4659763 * 0.3301256 *
*        4 *        hu * flowrateM * 21.192361 * 21.192361 * 8.8498989 *
*        5 *        hl * flowrateM * -32.74887 * 32.748874 * 30.146134 *
*        6 *         l * flowrateM * -23.89328 * 23.893280 * 8.2781934 *
*        7 *        kw * flowrateM * 7.5766167 * 7.5766167 * 2.7457665 *
************************************************************************