{
  "cells": [
    {
      "cell_type": "code",
      "execution_count": null,
      "metadata": {
        "collapsed": false
      },
      "outputs": [],
      "source": [
        "%matplotlib inline"
      ]
    },
    {
      "cell_type": "markdown",
      "metadata": {},
      "source": [
        "\n# Interface to PyKrige\n\nTo use fancier methods like\n[regression kriging](https://en.wikipedia.org/wiki/Regression-kriging)_,\nwe provide an interface to\n[PyKrige](https://github.com/GeoStat-Framework/PyKrige)_ (>v1.5), which means\nyou can pass a GSTools covariance model to the kriging routines of PyKrige.\n\nTo demonstrate the general workflow, we compare ordinary kriging of PyKrige\nwith the corresponding GSTools routine in 2D:\n"
      ]
    },
    {
      "cell_type": "code",
      "execution_count": null,
      "metadata": {
        "collapsed": false
      },
      "outputs": [],
      "source": [
        "import numpy as np\nfrom matplotlib import pyplot as plt\nfrom pykrige.ok import OrdinaryKriging\n\nimport gstools as gs\n\n# conditioning data\ncond_x = [0.3, 1.9, 1.1, 3.3, 4.7]\ncond_y = [1.2, 0.6, 3.2, 4.4, 3.8]\ncond_val = [0.47, 0.56, 0.74, 1.47, 1.74]\n\n# grid definition for output field\ngridx = np.arange(0.0, 5.5, 0.1)\ngridy = np.arange(0.0, 6.5, 0.1)"
      ]
    },
    {
      "cell_type": "markdown",
      "metadata": {},
      "source": [
        "A GSTools based :any:`Gaussian` covariance model:\n\n"
      ]
    },
    {
      "cell_type": "code",
      "execution_count": null,
      "metadata": {
        "collapsed": false
      },
      "outputs": [],
      "source": [
        "model = gs.Gaussian(\n    dim=2, len_scale=1, anis=0.2, angles=-0.5, var=0.5, nugget=0.1\n)"
      ]
    },
    {
      "cell_type": "markdown",
      "metadata": {},
      "source": [
        "## Ordinary Kriging with PyKrige\n\nOne can pass the defined GSTools model as\nvariogram model, which will `not` be fitted to the given data.\nBy providing the GSTools model, rotation and anisotropy are also\nautomatically defined:\n\n"
      ]
    },
    {
      "cell_type": "code",
      "execution_count": null,
      "metadata": {
        "collapsed": false
      },
      "outputs": [],
      "source": [
        "OK1 = OrdinaryKriging(cond_x, cond_y, cond_val, variogram_model=model)\nz1, ss1 = OK1.execute(\"grid\", gridx, gridy)\nplt.imshow(z1, origin=\"lower\")\nplt.show()"
      ]
    },
    {
      "cell_type": "markdown",
      "metadata": {},
      "source": [
        "## Ordinary Kriging with GSTools\n\nThe :any:`Ordinary` kriging class is provided by GSTools as a shortcut to\ndefine ordinary kriging with the general :any:`Krige` class.\n\nPyKrige's routines are using exact kriging by default (when given a nugget).\nTo reproduce this behavior in GSTools, we have to set ``exact=True``.\n\n"
      ]
    },
    {
      "cell_type": "code",
      "execution_count": null,
      "metadata": {
        "collapsed": false
      },
      "outputs": [],
      "source": [
        "OK2 = gs.krige.Ordinary(model, [cond_x, cond_y], cond_val, exact=True)\nOK2.structured([gridx, gridy])\nax = OK2.plot()\nax.set_aspect(\"equal\")"
      ]
    }
  ],
  "metadata": {
    "kernelspec": {
      "display_name": "Python 3",
      "language": "python",
      "name": "python3"
    },
    "language_info": {
      "codemirror_mode": {
        "name": "ipython",
        "version": 3
      },
      "file_extension": ".py",
      "mimetype": "text/x-python",
      "name": "python",
      "nbconvert_exporter": "python",
      "pygments_lexer": "ipython3",
      "version": "3.8.6"
    }
  },
  "nbformat": 4,
  "nbformat_minor": 0
}