.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "examples/10_normalizer/00_lognormal_kriging.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_examples_10_normalizer_00_lognormal_kriging.py: Log-Normal Kriging ------------------ Log Normal kriging is a term to describe a special workflow for kriging to deal with log-normal data, like conductivity or transmissivity in hydrogeology. It simply means to first convert the input data to a normal distribution, i.e. applying a logarithic function, then interpolating these values with kriging and transforming the result back with the exponential function. The resulting kriging variance describes the error variance of the log-values of the target variable. In this example we will use ordinary kriging. .. GENERATED FROM PYTHON SOURCE LINES 17-30 .. code-block:: Python import numpy as np import gstools as gs # condtions cond_pos = [0.3, 1.9, 1.1, 3.3, 4.7] cond_val = [0.47, 0.56, 0.74, 1.47, 1.74] # resulting grid gridx = np.linspace(0.0, 15.0, 151) # stable covariance model model = gs.Stable(dim=1, var=0.5, len_scale=2.56, alpha=1.9) .. GENERATED FROM PYTHON SOURCE LINES 31-33 In order to result in log-normal kriging, we will use the :any:`LogNormal` Normalizer. This is a parameter-less normalizer, so we don't have to fit it. .. GENERATED FROM PYTHON SOURCE LINES 33-35 .. code-block:: Python normalizer = gs.normalizer.LogNormal .. GENERATED FROM PYTHON SOURCE LINES 36-42 Now we generate the interpolated field as well as the mean field. This can be done by setting `only_mean=True` in :any:`Krige.__call__`. The result is then stored as `mean_field`. In terms of log-normal kriging, this mean represents the geometric mean of the field. .. GENERATED FROM PYTHON SOURCE LINES 42-48 .. code-block:: Python krige = gs.krige.Ordinary(model, cond_pos, cond_val, normalizer=normalizer) # interpolate the field krige(gridx) # also generate the mean field krige(gridx, only_mean=True) .. GENERATED FROM PYTHON SOURCE LINES 49-50 And that's it. Let's have a look at the results. .. GENERATED FROM PYTHON SOURCE LINES 50-56 .. code-block:: Python ax = krige.plot() # plotting the geometric mean krige.plot("mean_field", ax=ax) # plotting the conditioning data ax.scatter(cond_pos, cond_val, color="k", zorder=10, label="Conditions") ax.legend() .. image-sg:: /examples/10_normalizer/images/sphx_glr_00_lognormal_kriging_001.png :alt: Field 1D: (151,) :srcset: /examples/10_normalizer/images/sphx_glr_00_lognormal_kriging_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 0.054 seconds) .. _sphx_glr_download_examples_10_normalizer_00_lognormal_kriging.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: 00_lognormal_kriging.ipynb <00_lognormal_kriging.ipynb>` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: 00_lognormal_kriging.py <00_lognormal_kriging.py>` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: 00_lognormal_kriging.zip <00_lognormal_kriging.zip>` .. only:: html .. rst-class:: sphx-glr-signature `Gallery generated by Sphinx-Gallery `_