.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_calibration/bayesian_calibration/plot_bayesian_calibration_flooding.py" .. LINE NUMBERS ARE GIVEN BELOW. .. only:: html .. note:: :class: sphx-glr-download-link-note :ref:`Go to the end ` to download the full example code .. rst-class:: sphx-glr-example-title .. _sphx_glr_auto_calibration_bayesian_calibration_plot_bayesian_calibration_flooding.py: Bayesian calibration of the flooding model ========================================== .. GENERATED FROM PYTHON SOURCE LINES 6-10 Abstract -------- The goal of this example is to present the Bayesian calibration of the :ref:`flooding model`. .. GENERATED FROM PYTHON SOURCE LINES 13-74 Parameters to calibrate ----------------------- The vector of parameters to calibrate is: .. math:: \vect{\theta} = (K_s,Z_v,Z_m). The variables to calibrate are :math:`(K_s,Z_v,Z_m)` and are set to the following values: .. math:: K_s = 30, \qquad Z_v = 50, \qquad Z_m = 55. Observations ------------ In this section, we describe the statistical model associated with the :math:`n` observations. The errors of the water heights are associated with a gaussian distribution with a zero mean and a standard variation equal to: .. math:: \sigma=0.1. Therefore, the observed water heights are: .. math:: H_i = G(Q_i,K_s,Z_v,Z_m) + \epsilon_i for :math:`i=1,...,n` where .. math:: \epsilon \sim \mathcal{N}(0,\sigma^2) and we make the hypothesis that the observation errors are independent. We consider a sample size equal to: .. math:: n=20. The observations are the couples :math:`\{(Q_i,H_i)\}_{i=1,...,n}`, i.e. each observation is a couple made of the flowrate and the corresponding river height. Variables --------- - :math:`Q` : Input. Observed. - :math:`K_s`, :math:`Z_v`, :math:`Z_m` : Input. Calibrated. - :math:`H`: Output. Observed. Analysis -------- In the description of the :ref:`flooding model`, we see that only one parameter can be identified. Hence, calibrating this model requires some regularization. In this example, we use Bayesian methods as a way to regularize the model. .. GENERATED FROM PYTHON SOURCE LINES 76-78 Generate the observations ------------------------- .. GENERATED FROM PYTHON SOURCE LINES 80-88 .. code-block:: Python import pylab as pl from openturns.usecases import flood_model import openturns.viewer as viewer import numpy as np import openturns as ot ot.Log.Show(ot.Log.NONE) .. GENERATED FROM PYTHON SOURCE LINES 89-90 A basic implementation of the probabilistic model is available in the usecases module : .. GENERATED FROM PYTHON SOURCE LINES 90-92 .. code-block:: Python fm = flood_model.FloodModel() .. GENERATED FROM PYTHON SOURCE LINES 93-101 We define the model :math:`g` which has 4 inputs and one output H. The nonlinear least squares does not take into account for bounds in the parameters. Therefore, we ensure that the output is computed whatever the inputs. The model fails into two situations: * if :math:`K_s<0`, * if :math:`Z_v-Z_m<0`. In these cases, we return an infinite number. .. GENERATED FROM PYTHON SOURCE LINES 103-117 .. code-block:: Python def functionFlooding(X): L = 5.0e3 B = 300.0 Q, K_s, Z_v, Z_m = X alpha = (Z_m - Z_v) / L if alpha < 0.0 or K_s <= 0.0: H = np.inf else: H = (Q / (K_s * B * np.sqrt(alpha))) ** (3.0 / 5.0) return [H] .. GENERATED FROM PYTHON SOURCE LINES 118-122 .. code-block:: Python g = ot.PythonFunction(4, 1, functionFlooding) g = ot.MemoizeFunction(g) g.setOutputDescription(["H (m)"]) .. GENERATED FROM PYTHON SOURCE LINES 123-124 We load the input distribution for :math:`Q`. .. GENERATED FROM PYTHON SOURCE LINES 124-126 .. code-block:: Python Q = fm.Q .. GENERATED FROM PYTHON SOURCE LINES 127-128 Set the parameters to be calibrated. .. GENERATED FROM PYTHON SOURCE LINES 130-137 .. code-block:: Python K_s = ot.Dirac(30.0) Z_v = ot.Dirac(50.0) Z_m = ot.Dirac(55.0) K_s.setDescription(["Ks (m^(1/3)/s)"]) Z_v.setDescription(["Zv (m)"]) Z_m.setDescription(["Zm (m)"]) .. GENERATED FROM PYTHON SOURCE LINES 138-139 We create the joint input distribution. .. GENERATED FROM PYTHON SOURCE LINES 141-143 .. code-block:: Python inputRandomVector = ot.ComposedDistribution([Q, K_s, Z_v, Z_m]) .. GENERATED FROM PYTHON SOURCE LINES 144-145 Create a Monte-Carlo sample of the output :math:`H`. .. GENERATED FROM PYTHON SOURCE LINES 147-151 .. code-block:: Python nbobs = 20 inputSample = inputRandomVector.getSample(nbobs) outputH = g(inputSample) .. GENERATED FROM PYTHON SOURCE LINES 152-153 Generate the observation noise and add it to the output of the model. .. GENERATED FROM PYTHON SOURCE LINES 155-161 .. code-block:: Python sigmaObservationNoiseH = 0.1 # (m) noiseH = ot.Normal(0.0, sigmaObservationNoiseH) ot.RandomGenerator.SetSeed(0) sampleNoiseH = noiseH.getSample(nbobs) Hobs = outputH + sampleNoiseH .. GENERATED FROM PYTHON SOURCE LINES 162-163 Plot the Y observations versus the X observations. .. GENERATED FROM PYTHON SOURCE LINES 165-167 .. code-block:: Python Qobs = inputSample[:, 0] .. GENERATED FROM PYTHON SOURCE LINES 168-174 .. code-block:: Python graph = ot.Graph("Observations", "Q (m3/s)", "H (m)", True) cloud = ot.Cloud(Qobs, Hobs) graph.add(cloud) view = viewer.View(graph) .. image-sg:: /auto_calibration/bayesian_calibration/images/sphx_glr_plot_bayesian_calibration_flooding_001.png :alt: Observations :srcset: /auto_calibration/bayesian_calibration/images/sphx_glr_plot_bayesian_calibration_flooding_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 175-177 Setting the calibration parameters ---------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 179-184 Define the parametric model :math:`\vect z = f_Q(\vect\theta)` that associates each observation :math:`Q` and value of the parameters :math:`\vect \theta = (K_s, Z_v, Z_m)` to the parameters of the distribution of the corresponding observation: here :math:`\vect z=(\mu, \sigma)` with :math:`\mu = G(Q, K_s, Z_v, Z_m)` and :math:`\sigma = 0.5`. .. GENERATED FROM PYTHON SOURCE LINES 186-197 .. code-block:: Python def fullModelPy(X): Q, K_s, Z_v, Z_m = X mu = g(X)[0] sigma = 0.5 # (m^2) The standard deviation of the observation error. return [mu, sigma] fullModel = ot.PythonFunction(4, 2, fullModelPy) linkFunction = ot.ParametricFunction(fullModel, [0], [np.nan]) print(linkFunction) .. rst-class:: sphx-glr-script-out .. code-block:: none ParametricEvaluation(class=PythonEvaluation name=OpenTURNSPythonFunction, parameters positions=[0], parameters=[x0 : nan], input positions=[1,2,3]) .. GENERATED FROM PYTHON SOURCE LINES 198-200 Define the value of the reference values of the :math:`\vect\theta` parameter. In the Bayesian framework, this is called the mean of the *prior* Gaussian distribution. In the data assimilation framework, this is called the *background*. .. GENERATED FROM PYTHON SOURCE LINES 202-208 .. code-block:: Python KsInitial = 20.0 ZvInitial = 49.0 ZmInitial = 51.0 parameterPriorMean = [KsInitial, ZvInitial, ZmInitial] paramDim = len(parameterPriorMean) .. GENERATED FROM PYTHON SOURCE LINES 209-210 Define the covariance matrix of the parameters :math:`\vect\theta` to calibrate. .. GENERATED FROM PYTHON SOURCE LINES 212-216 .. code-block:: Python sigmaKs = 5.0 sigmaZv = 1.0 sigmaZm = 1.0 .. GENERATED FROM PYTHON SOURCE LINES 217-222 .. code-block:: Python parameterPriorCovariance = ot.CovarianceMatrix(paramDim) parameterPriorCovariance[0, 0] = sigmaKs ** 2 parameterPriorCovariance[1, 1] = sigmaZv ** 2 parameterPriorCovariance[2, 2] = sigmaZm ** 2 .. GENERATED FROM PYTHON SOURCE LINES 223-224 Define the prior distribution :math:`\pi(\vect\theta)` of the parameter :math:`\vect\theta` .. GENERATED FROM PYTHON SOURCE LINES 226-229 .. code-block:: Python prior = ot.Normal(parameterPriorMean, parameterPriorCovariance) prior.setDescription(["Ks", "Zv", "Zm"]) .. GENERATED FROM PYTHON SOURCE LINES 230-235 Define the distribution of observations :math:`\vect{y} | \vect{z}` conditional on model predictions. Note that its parameter dimension is the one of :math:`\vect{z}`, so the model must be adjusted accordingly. In other words, the input argument of the `setParameter` method of the conditional distribution must be equal to the dimension of the output of the `model`. Hence, we do not have to set the actual parameters: only the type of distribution is used. .. GENERATED FROM PYTHON SOURCE LINES 237-239 .. code-block:: Python conditional = ot.Normal() .. GENERATED FROM PYTHON SOURCE LINES 240-243 The proposed steps for :math:`K_s` :math:`Z_v` and :math:`Z_m` will all follow uniform distributions, but with different supports. .. GENERATED FROM PYTHON SOURCE LINES 245-247 .. code-block:: Python proposal = [ot.Uniform(-5.0, 5.0), ot.Uniform(-1.0, 1.0), ot.Uniform(-1.0, 1.0)] .. GENERATED FROM PYTHON SOURCE LINES 248-250 Build a Gibbs sampler --------------------- .. GENERATED FROM PYTHON SOURCE LINES 250-260 .. code-block:: Python initialState = parameterPriorMean mh_coll = [ ot.RandomWalkMetropolisHastings(prior, initialState, proposal[i], [i]) for i in range(paramDim) ] for mh in mh_coll: mh.setLikelihood(conditional, Hobs, linkFunction, Qobs) sampler = ot.Gibbs(mh_coll) .. GENERATED FROM PYTHON SOURCE LINES 261-262 Generate a sample from the posterior distribution of :math:`\vect \theta`. .. GENERATED FROM PYTHON SOURCE LINES 264-267 .. code-block:: Python sampleSize = 1000 sample = sampler.getSample(sampleSize) .. GENERATED FROM PYTHON SOURCE LINES 268-269 Look at the acceptance rates of the random walk Metropolis-Hastings samplers. .. GENERATED FROM PYTHON SOURCE LINES 269-272 .. code-block:: Python [mh.getAcceptanceRate() for mh in sampler.getMetropolisHastingsCollection()] .. rst-class:: sphx-glr-script-out .. code-block:: none [0.365, 0.359, 0.333] .. GENERATED FROM PYTHON SOURCE LINES 273-274 Build the distribution of the posterior by kernel smoothing. .. GENERATED FROM PYTHON SOURCE LINES 274-278 .. code-block:: Python kernel = ot.KernelSmoothing() posterior = kernel.build(sample) .. GENERATED FROM PYTHON SOURCE LINES 279-280 Display prior vs posterior for each parameter. .. GENERATED FROM PYTHON SOURCE LINES 280-293 .. code-block:: Python fig = pl.figure(figsize=(12, 4)) for parameter_index in range(paramDim): graph = posterior.getMarginal(parameter_index).drawPDF() priorGraph = prior.getMarginal(parameter_index).drawPDF() priorGraph.setColors(["blue"]) graph.add(priorGraph) graph.setLegends(["Posterior", "Prior"]) ax = fig.add_subplot(1, paramDim, parameter_index + 1) _ = ot.viewer.View(graph, figure=fig, axes=[ax]) _ = fig.suptitle("Bayesian calibration") .. image-sg:: /auto_calibration/bayesian_calibration/images/sphx_glr_plot_bayesian_calibration_flooding_002.png :alt: Bayesian calibration :srcset: /auto_calibration/bayesian_calibration/images/sphx_glr_plot_bayesian_calibration_flooding_002.png :class: sphx-glr-single-img .. _sphx_glr_download_auto_calibration_bayesian_calibration_plot_bayesian_calibration_flooding.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_bayesian_calibration_flooding.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_bayesian_calibration_flooding.py `