File: plot_kriging_multioutput_firesatellite.py

package info (click to toggle)
openturns 1.24-4
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 66,204 kB
  • sloc: cpp: 256,662; python: 63,381; ansic: 4,414; javascript: 406; sh: 180; xml: 164; yacc: 123; makefile: 98; lex: 55
file content (131 lines) | stat: -rw-r--r-- 4,061 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
"""
Example of multi output Kriging on the fire satellite model
===========================================================
"""

# %%
# This example aims to illustrate Kriging metamodel with several outputs on the fire satellite model.


# %%
# Loading of the model
# --------------------
# This model involves 9 input variables and 3 output variables.
# We load the :ref:`Fire satellite use case<use-case-firesatellite>`.

# %%
import openturns as ot
from openturns.usecases.fireSatellite_function import FireSatelliteModel
from openturns.viewer import View

ot.Log.Show(ot.Log.NONE)

m = FireSatelliteModel()

# %%
# We define the function that evaluates the outputs depending on the inputs.

# %%
model = m.model

# %%
# We also define the distribution of input variables to build the training and test sets.

# %%
inputDistribution = m.inputDistribution


# %%
# Generation of data
# ------------------
# We now generate the input and output training sets as 10 times the dimension of the input vector.

# %%
ot.RandomGenerator.SetSeed(0)
experiment = ot.LHSExperiment(inputDistribution, 10 * m.dim)
inputTrainingSet = experiment.generate()
outputTrainingSet = model(inputTrainingSet)

print("Lower and upper bounds of inputTrainingSet:")
print(inputTrainingSet.getMin(), inputTrainingSet.getMax())

# %%
# Creation of metamodel
# ---------------------
# We choose to use a constant trend.

# %%
linear_basis = ot.LinearBasisFactory(m.dim).build()
basis = ot.Basis(
    [
        ot.AggregatedFunction([linear_basis.build(k)] * 3)
        for k in range(linear_basis.getSize())
    ]
)

# %%
# We would like to have separate covariance models for the three outputs.
# To do so, we use the :class:`~openturns.TensorizedCovarianceModel`.
# For the purpose of illustration, we consider :class:`~openturns.MaternModel` for the first and third outputs, and :class:`~openturns.SquaredExponential` for the second output.

# %%
myCov1 = ot.MaternModel([1.0] * m.dim, 2.5)
myCov2 = ot.SquaredExponential([1.0] * m.dim)
myCov3 = ot.MaternModel([1.0] * m.dim, 2.5)

covarianceModel = ot.TensorizedCovarianceModel([myCov1, myCov2, myCov3])

# %%
# The scaling of the data is really important when dealing with Kriging,
# especially considering the domain definition of the input variables (the
# altitude varies in order of :math:`10^7` whereas the drag coefficient is around 1).
# We thus define appropriate bounds for the training algorithm based on the
# domain definition of each variable.

# %%
scaleOptimizationBounds = ot.Interval(
    [1.0e6, 1.0e3, 1.0e3, 1.0, 0.1, 0.1, 0.1, 0.1, 0.1, 0.1, 0.1, 0.1],
    [2.0e7, 2.0e3, 2.0e3, 1e2, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0, 10.0],
)

# %%
# We can now define the scaled version of Kriging model.
algo = ot.KrigingAlgorithm(inputTrainingSet, outputTrainingSet, covarianceModel, basis)
algo.setOptimizationBounds(scaleOptimizationBounds)
algo.setOptimizeParameters(True)

# %%
# We run the algorithm and get the metamodel.
algo.run()
result = algo.getResult()
krigingMetamodel = result.getMetaModel()

# %%
# Validation of metamodel
# -----------------------
# To validate the metamodel, we create a validation set of size equal to 50 times the input vector dimension to evaluate the functions.

# %%
ot.RandomGenerator.SetSeed(1)
experimentTest = ot.LHSExperiment(inputDistribution, 50 * m.dim)
inputTestSet = experimentTest.generate()
outputTestSet = model(inputTestSet)

# %%
# Then, we use the :class:`~openturns.MetaModelValidation` class to validate the metamodel.
metamodelPredictions = krigingMetamodel(inputTestSet)
val = ot.MetaModelValidation(outputTestSet, metamodelPredictions)

r2Score = val.computeR2Score()

label = ["Total torque", "Total power", "Solar array area"]

for i in range(3):
    graph = val.drawValidation().getGraph(0, i)
    graph.setLegends([""])
    graph.setLegends(["R2 = %.2f%%" % (100 * r2Score[i]), ""])
    graph.setLegendPosition("upper left")
    graph.setXTitle("Exact function")
    graph.setYTitle("Metamodel prediction")
    graph.setTitle(label[i])
    View(graph)