13.4.3. Macro “samplingLHSCorrelation.py

13.4.3.1. Objective

Generate a design-of-experiments, of 3000 patterns using the LHS method, with 3 random attributes and taking into account a correlation between the first two attributes. This correlation, equals to 0.99, is computed using the rank values (the Spearmann definition, c.f. Description of a correlation and in [Bla17] for more details).

  1. Attribute x1 follows an uniform law on interval [3, 4];

  2. Attribute x2 follows a normal law with a mean value equal to 0.5 and 1.5 as standard deviation;

  3. Attribute x3 obeys a triangular law on interval [1, 5] with mode 4.

13.4.3.2. Macro Uranie

"""
Example of simple LHS correlated DoE generation
"""
from URANIE import DataServer, Sampler
import ROOT

# Create a DataServer.TDataServer
tds = DataServer.TDataServer()
# Fill the DataServer with the three attributes of the study
tds.addAttribute(DataServer.TUniformDistribution("x1", 3., 4.))
tds.addAttribute(DataServer.TNormalDistribution("x2", 0.5, 1.5))
tds.addAttribute(DataServer.TTriangularDistribution("x3", 1., 5., 4.))

# Generate the sampling from the DataServer.TDataServer
sampling = Sampler.TSampling(tds, "lhs", 3000)
sampling.setUserCorrelation(0, 1, 0.99)
sampling.generateSample()

tds.exportData("toto.dat", "x1:x3")
tds.exportDataHeader("toto.h", "x1:x2")

# Graph
Canvas = ROOT.TCanvas("c1", "Graph samplingCorrelation", 5, 64, 1270, 667)
pad = ROOT.TPad("pad", "pad", 0, 0.03, 1, 1)
pad.Draw()
pad.Divide(2, 2)

pad.cd(1)
tds.Draw("x1")
pad.cd(2)
tds.Draw("x2")
pad.cd(3)
tds.Draw("x3")
pad.cd(4)
tds.drawTufte("x1:x2")

Laws are set for each attribute and linked to a TDataServer:

tds.addAttribute(DataServer.TUniformDistribution("x1", 3., 4.))
tds.addAttribute(DataServer.TNormalDistribution("x2", 0.5, 1.5))
tds.addAttribute(DataServer.TTriangularDistribution("x3", 1., 5., 4.))

The sampling is generated with the LHS method and a correlation is set between the two first attributes:

sampling = Sampler.TSampling(tds, "lhs", 3000)
sampling.setUserCorrelation(0, 1, 0.99)
sampling.generateSample()

Data are partially exported in an ASCII file and a header file:

tds.exportData("toto.dat", "x1:x3")
tds.exportDataHeader("toto.h", "x1:x2")

13.4.3.3. Graph

../../_images/samplingLHSCorrelation.png

Figure 13.11 Graph de la macro “samplingLHSCorrelation.py”