{ "cells": [ { "cell_type": "markdown", "metadata": {}, "source": [ "\n# Elder\n\nThe classic 2D Elder problem demonstrates free convection.\nTraditionally this was created for heat transport, but we use\na modified version for salt transport.\nThe conceptual model can be seen as a 2D sand box,\nwith on top a salt lake in the center and fresh lakes\non both the outer edges of the top row.\n\nMore info about the theory behind the Elder problem:\n\nSimpson, J., & Clement, P. (2003).\nTheoretical analysis of the worthiness of Henry and Elder\nproblems as benchmark of density-dependent groundwater flow models.\n`Advances in Water Resources, 1708` (02).\nRetrieved from http://www.eng.auburn.edu/~clemept/publsihed_pdf/awrmat.pdf\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 matplotlib.pyplot as plt\nimport numpy as np\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 = 160\nnlay = 82\n\ndz = 1.875\ndx = 3.75\ndy = -dx" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "setup tops and bottoms\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "top1D = xr.DataArray(\n np.arange(nlay * dz, 0.0, -dz), {\"layer\": np.arange(1, nlay + 1)}, (\"layer\")\n)\n\nbot = top1D - dz\ntop = nlay * dz" ] }, { "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\nfig, ax = plt.subplots()\nbnd.plot(y=\"layer\", yincrease=False, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Boundary Conditions\n\nSet the constant heads by specifying a negative value in iboud,\nthat is: ``bnd[index] = -1```\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "bnd[0, :, 0:40] = 0\nbnd[0, :, 121:160] = 0\nbnd[1, :, 0] = -1\nbnd[1, :, 159] = -1\n\nfig, ax = plt.subplots()\nbnd.plot(y=\"layer\", yincrease=False, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "Define the icbund, which sets which cells\nin the solute transport model are active, inactive or constant.\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "icbund = 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, nlay + 1),\n },\n dims=(\"layer\", \"y\", \"x\"),\n)\n\nicbund[81, :, :] = -1\nicbund[0, :, 41:120] = -1\n\nfig, ax = plt.subplots()\nicbund.plot(y=\"layer\", yincrease=False, ax=ax)" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Initial Conditions\n\nDefine the starting concentration\n\n" ] }, { "cell_type": "code", "execution_count": null, "metadata": { "collapsed": false }, "outputs": [], "source": [ "sconc = xr.DataArray(\n data=np.full((nlay, nrow, ncol), 0.0),\n coords={\n \"y\": [0.5],\n \"x\": np.arange(0.5 * dx, dx * ncol, dx),\n \"layer\": np.arange(1, nlay + 1),\n },\n dims=(\"layer\", \"y\", \"x\"),\n)\n\nsconc[81, :, :] = 0\nsconc[0, :, 41:120] = 280.0\n\nfig, ax = plt.subplots()\nsconc.plot(y=\"layer\", yincrease=False, ax=ax)" ] }, { "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(\"Elder\")\nm[\"bas\"] = imod.wq.BasicFlow(ibound=bnd, top=top, bottom=bot, starting_head=0.0)\nm[\"lpf\"] = imod.wq.LayerPropertyFlow(\n k_horizontal=0.411, k_vertical=0.411, specific_storage=0.0\n)\nm[\"btn\"] = imod.wq.BasicTransport(\n icbund=icbund, starting_concentration=sconc, porosity=0.1\n)\nm[\"adv\"] = imod.wq.AdvectionTVD(courant=1.0)\nm[\"dsp\"] = imod.wq.Dispersion(longitudinal=0.0, diffusion_coefficient=0.308)\nm[\"vdf\"] = imod.wq.VariableDensityFlow(density_concentration_slope=0.71)\nm[\"wel\"] = imod.wq.Well(id_name=\"wel\", x=0.5 * dx, y=0.5, rate=0.28512)\nm[\"pcg\"] = imod.wq.PreconditionedConjugateGradientSolver(\n max_iter=150, inner_iter=30, hclose=0.0001, rclose=0.1, 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(additional_times=[\"2000-01-01T00:00\", \"2020-01-01T00:00\"])" ] }, { "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```batch\ncd c:\\path\\to\\modeldir\nc:\\path\\to\\imod\\folder\\iMOD-WQ_V5_3_SVN359_X64R.exe Elder.run\n```\nNote that the version name of your executable might differ.\n\n\n" ] }, { "cell_type": "markdown", "metadata": {}, "source": [ "## Visualise results\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.11.6" } }, "nbformat": 4, "nbformat_minor": 0 }