.. 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 Click :ref:`here ` 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-87 .. code-block:: default 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 88-89 A basic implementation of the probabilistic model is available in the usecases module : .. GENERATED FROM PYTHON SOURCE LINES 89-91 .. code-block:: default fm = flood_model.FloodModel() .. GENERATED FROM PYTHON SOURCE LINES 92-100 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 102-116 .. code-block:: default 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 117-121 .. code-block:: default g = ot.PythonFunction(4, 1, functionFlooding) g = ot.MemoizeFunction(g) g.setOutputDescription(["H (m)"]) .. GENERATED FROM PYTHON SOURCE LINES 122-123 We load the input distribution for :math:`Q`. .. GENERATED FROM PYTHON SOURCE LINES 123-125 .. code-block:: default Q = fm.Q .. GENERATED FROM PYTHON SOURCE LINES 126-127 Set the parameters to be calibrated. .. GENERATED FROM PYTHON SOURCE LINES 129-136 .. code-block:: default 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 137-138 We create the joint input distribution. .. GENERATED FROM PYTHON SOURCE LINES 140-142 .. code-block:: default inputRandomVector = ot.ComposedDistribution([Q, K_s, Z_v, Z_m]) .. GENERATED FROM PYTHON SOURCE LINES 143-144 Create a Monte-Carlo sample of the output :math:`H`. .. GENERATED FROM PYTHON SOURCE LINES 146-150 .. code-block:: default nbobs = 20 inputSample = inputRandomVector.getSample(nbobs) outputH = g(inputSample) .. GENERATED FROM PYTHON SOURCE LINES 151-152 Generate the observation noise and add it to the output of the model. .. GENERATED FROM PYTHON SOURCE LINES 154-160 .. code-block:: default sigmaObservationNoiseH = 0.1 # (m) noiseH = ot.Normal(0., sigmaObservationNoiseH) ot.RandomGenerator.SetSeed(0) sampleNoiseH = noiseH.getSample(nbobs) Hobs = outputH + sampleNoiseH .. GENERATED FROM PYTHON SOURCE LINES 161-162 Plot the Y observations versus the X observations. .. GENERATED FROM PYTHON SOURCE LINES 164-166 .. code-block:: default Qobs = inputSample[:, 0] .. GENERATED FROM PYTHON SOURCE LINES 167-173 .. code-block:: default 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 174-176 Setting the calibration parameters ---------------------------------- .. GENERATED FROM PYTHON SOURCE LINES 178-183 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 185-196 .. code-block:: default 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 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 197-198 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 200-206 .. code-block:: default KsInitial = 20. ZvInitial = 49. ZmInitial = 51. parameterPriorMean = [KsInitial, ZvInitial, ZmInitial] paramDim = len(parameterPriorMean) .. GENERATED FROM PYTHON SOURCE LINES 207-208 Define the covariance matrix of the parameters :math:`\vect\theta` to calibrate. .. GENERATED FROM PYTHON SOURCE LINES 210-214 .. code-block:: default sigmaKs = 5. sigmaZv = 1. sigmaZm = 1. .. GENERATED FROM PYTHON SOURCE LINES 215-220 .. code-block:: default parameterPriorCovariance = ot.CovarianceMatrix(paramDim) parameterPriorCovariance[0, 0] = sigmaKs**2 parameterPriorCovariance[1, 1] = sigmaZv**2 parameterPriorCovariance[2, 2] = sigmaZm**2 .. GENERATED FROM PYTHON SOURCE LINES 221-222 Define the the prior distribution :math:`\pi(\vect\theta)` of the parameter :math:`\vect\theta` .. GENERATED FROM PYTHON SOURCE LINES 224-227 .. code-block:: default prior = ot.Normal(parameterPriorMean, parameterPriorCovariance) prior.setDescription(['Ks', 'Zv', 'Zm']) .. GENERATED FROM PYTHON SOURCE LINES 228-231 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 233-235 .. code-block:: default conditional = ot.Normal() .. GENERATED FROM PYTHON SOURCE LINES 236-239 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 241-243 .. code-block:: default proposal = [ot.Uniform(-5., 5.), ot.Uniform(-1., 1.), ot.Uniform(-1., 1.)] .. GENERATED FROM PYTHON SOURCE LINES 244-246 Build a Gibbs sampler --------------------- .. GENERATED FROM PYTHON SOURCE LINES 246-254 .. code-block:: default 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 255-256 Tuning of the Gibbs algorithm. .. GENERATED FROM PYTHON SOURCE LINES 258-261 .. code-block:: default sampler.setThinning(1) sampler.setBurnIn(200) .. GENERATED FROM PYTHON SOURCE LINES 262-263 Generate a sample from the posterior distribution of :math:`\vect \theta`. .. GENERATED FROM PYTHON SOURCE LINES 265-268 .. code-block:: default sampleSize = 1000 sample = sampler.getSample(sampleSize) .. GENERATED FROM PYTHON SOURCE LINES 269-270 Look at the acceptance rates of the random walk Metropolis-Hastings samplers. .. GENERATED FROM PYTHON SOURCE LINES 270-273 .. code-block:: default [mh.getAcceptanceRate() for mh in sampler.getMetropolisHastingsCollection()] .. rst-class:: sphx-glr-script-out Out: .. code-block:: none [0.5433333333333333, 0.655, 0.6416666666666667] .. GENERATED FROM PYTHON SOURCE LINES 274-275 Build the distribution of the posterior by kernel smoothing. .. GENERATED FROM PYTHON SOURCE LINES 275-279 .. code-block:: default kernel = ot.KernelSmoothing() posterior = kernel.build(sample) .. GENERATED FROM PYTHON SOURCE LINES 280-281 Display prior vs posterior for each parameter. .. GENERATED FROM PYTHON SOURCE LINES 281-294 .. code-block:: default 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 .. rst-class:: sphx-glr-timing **Total running time of the script:** ( 0 minutes 1.436 seconds) .. _sphx_glr_download_auto_calibration_bayesian_calibration_plot_bayesian_calibration_flooding.py: .. only :: html .. container:: sphx-glr-footer :class: sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_bayesian_calibration_flooding.py ` .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_bayesian_calibration_flooding.ipynb ` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_