{ "cells": [ { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "%matplotlib inline" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Henry\n\nThe classic 2D Henry problem demonstrates the development of a fresh-salt\ninterface.\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "We'll start with the usual imports\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "import numpy as np\nimport pandas as pd\nimport xarray as xr\n\nimport imod" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Discretization\n\nWe'll start off by creating a model discretization, since\nthis is a simple conceptual model.\nThe model is a 2D cross-section, hence ``nrow = 1``.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "nrow = 1\nncol = 100\nnlay = 50\n\ndz = 1.0\ndx = 1.0\ndy = -dx\n\ntop1D = xr.DataArray(\n np.arange(nlay * dz, 0.0, -dz), {\"layer\": np.arange(1, 1 + nlay)}, (\"layer\")\n)\nbot = top1D - 1.0" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Set up ibound, which sets where active cells are `(ibound = 1.0)`\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "bnd = xr.DataArray(\n data=np.full((nlay, nrow, ncol), 1.0),\n coords={\n \"y\": [0.5],\n \"x\": np.arange(0.5 * dx, dx * ncol, dx),\n \"layer\": np.arange(1, 1 + nlay),\n \"dx\": dx,\n \"dy\": dy,\n },\n dims=(\"layer\", \"y\", \"x\"),\n)\n\n# Boundary Conditions\n# -------------------\n#\n# We define constant head here, after generating the tops, or we'd end up with negative top values\nbnd[:, :, -1] = -1" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Create WEL package\n\nFirst we scale the pumping rate with discretization\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "qscaled = 0.03 * (dz * abs(dy))" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Fresh water injection with well\nAdd the arguments as a list, so pandas doesn't complain about having to set\nan index.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "weldata = pd.DataFrame()\nweldata[\"x\"] = [0.5]\nweldata[\"y\"] = [0.5]\nweldata[\"q\"] = [qscaled]" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Build\n\nFinally, we build the model.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "m = imod.wq.SeawatModel(\"Henry\")\nm[\"bas\"] = imod.wq.BasicFlow(ibound=bnd, top=50.0, bottom=bot, starting_head=1.0)\nm[\"lpf\"] = imod.wq.LayerPropertyFlow(\n k_horizontal=10.0, k_vertical=10.0, specific_storage=0.0\n)\nm[\"btn\"] = imod.wq.BasicTransport(\n icbund=bnd, starting_concentration=35.0, porosity=0.35\n)\nm[\"adv\"] = imod.wq.AdvectionTVD(courant=1.0)\nm[\"dsp\"] = imod.wq.Dispersion(longitudinal=0.1, diffusion_coefficient=1.0e-9)\nm[\"vdf\"] = imod.wq.VariableDensityFlow(density_concentration_slope=0.71)\nm[\"wel\"] = imod.wq.Well(\n id_name=\"well\", x=weldata[\"x\"], y=weldata[\"y\"], rate=weldata[\"q\"]\n)\nm[\"pcg\"] = imod.wq.PreconditionedConjugateGradientSolver(\n max_iter=150, inner_iter=30, hclose=0.0001, rclose=1.0, relax=0.98, damp=1.0\n)\nm[\"gcg\"] = imod.wq.GeneralizedConjugateGradientSolver(\n max_iter=150,\n inner_iter=30,\n cclose=1.0e-6,\n preconditioner=\"mic\",\n lump_dispersion=True,\n)\nm[\"oc\"] = imod.wq.OutputControl(save_head_idf=True, save_concentration_idf=True)\nm.create_time_discretization(\n additional_times=pd.date_range(\"2000-01-01\", \"2001-01-01\", freq=\"M\")\n)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Now we write the model\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "modeldir = imod.util.temporary_directory()\nm.write(modeldir, resultdir_is_workdir=True)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Run\n\nYou can run the model using the comand prompt and the iMOD-WQ executable.\nThis is part of the iMOD v5 release, which can be downloaded here:\nhttps://oss.deltares.nl/web/imod/download-imod5 .\nThis only works on Windows.\n\nTo run your model, open up a command prompt\nand run the following commands:\n\n.. code-block:: batch\n\n cd c:\\path\\to\\modeldir\n c:\\path\\to\\imod\\folder\\iMOD-WQ_V5_3_SVN359_X64R.exe Henry.run\n\nNote that the version name of your executable might differ.\n\n%%\nVisualise results\n-----------------\n\nAfter succesfully running the model, you can\nplot results as follows:\n\n.. code:: python\n\n head = imod.idf.open(modeldir / \"results/head/*.idf\")\n\n fig, ax = plt.subplots()\n head.plot(yincrease=False, ax=ax)\n\n conc = imod.idf.open(modeldir / \"results/conc/*.idf\")\n\n fig, ax = plt.subplots()\n conc.plot(levels=range(0, 35, 5), yincrease=False, ax=ax)\n\n\n" ] } ], "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.10.5" } }, "nbformat": 4, "nbformat_minor": 0 }