{ "cells": [ { "cell_type": "markdown", "id": "aaa60d19", "metadata": {}, "source": [ "# Reactive transport modeling along a rock core after injection of the fluid-rock composition\n", "\n", "In this tutorial, we show how Reaktoro can be used for sequential calculations of the reactive transport along a\n", "rock column after injecting the fluid and rock composition at temperature 60 °C and pressure 100 bar.\n", "\n", "Using Reaktoro in Python requires first an import of the python package **reaktoro**. From this point on,\n", "we can use the library components of Reaktoro (classes, methods, constants), which are needed to define our\n", "chemical system and chemical reaction modeling problems.\n", "\n", "## Importing python packages\n", "\n", "First, we need to import a few Python packages to enable us to perform the numerical calculations and plotting." ] }, { "cell_type": "code", "execution_count": null, "id": "725397a4", "metadata": {}, "outputs": [], "source": [ "print('============================================================')\n", "print('Make sure you have the following Python packages installed: ')\n", "print(' numpy, natsort, bokeh')\n", "print('These can be installed with pip:')\n", "print(' pip install numpy natsort bokeh')\n", "print('============================================================')\n", "from reaktoro import *\n", "import numpy as np\n", "from natsort import natsorted\n", "from tqdm.notebook import tqdm\n", "import os\n", "\n", "# Import components of bokeh library\n", "from bokeh.io import show, output_notebook\n", "from bokeh.layouts import column\n", "from bokeh.plotting import figure\n", "from bokeh.models import Range1d, ColumnDataSource\n", "from bokeh.layouts import gridplot" ] }, { "cell_type": "markdown", "id": "37d24b70", "metadata": {}, "source": [ "We import the **reaktoro** Python package so that we can use its classes and methods for performing chemical\n", "reaction calculations, **numpy** for working with arrays, **tqdm** for the progress bar functionality and **os**,\n", "to provide a portable way of using operating system dependent functionality. For plotting capabilities of obtained\n", "results, we use **bokeh** library.\n", "\n", "> **Note**: To make sure that all the widgets are working correctly, make sure to run:\n", "> `$ jupyter nbextension enable --py widgetsnbextension` and\n", "> `$jupyter labextension install @jupyter-widgets/jupyterlab-manager`.\n", "\n", "> **Note**: To simplify the tutorials, we use `from reaktoro import *`, which imports all components of the\n", "> **reaktoro** package into the default Python namespace. We note that this can potentially create name conflicts\n", "> when used in bigger projects. For your applications, consider using `import reaktoro as rkt` instead,\n", "> and call classes and methods as `rkt.Database`, `rkt.ChemicalSystem`, `rkt.equilibrate`, etc.\n", "\n", "## Initializing auxiliary time-related constants\n", "In this step, we initialize auxiliary time-related constants from seconds up to years used in the rest of the code." ] }, { "cell_type": "code", "execution_count": null, "id": "bd41e1e7", "metadata": {}, "outputs": [], "source": [ "second = 1\n", "minute = 60\n", "hour = 60 * minute\n", "day = 24 * hour\n", "week = 7 * day\n", "year = 365 * day" ] }, { "cell_type": "markdown", "id": "08c1b802", "metadata": {}, "source": [ "## Defining parameters for the reactive transport simulation\n", "\n", "Next, we define reactive transport and numerical discretization parameters. In particular, we specify the considered\n", "rock domain by setting coordinates of its left and right boundaries to 0.0 m and 100.0 m, respectively. The\n", "discretization parameters, i.e., the number of cells and steps in time, are both set to 100. The reactive\n", "transport modeling procedure assumes a constant fluid velocity of 1 m/week (1.65 · 10-6 m/s) and\n", "the same diffusion coefficient of 10-9 m2/s for all fluid species (without dispersivity).\n", "The size of the time-step is set to 30 minutes. Temperature and pressure are set to 60 °C and 100 bar,\n", "respectively, throughout the whole tutorial." ] }, { "cell_type": "code", "execution_count": null, "id": "ec8db0ef", "metadata": {}, "outputs": [], "source": [ "# Discretization parameters\n", "xl = 0.0 # x-coordinate of the left boundary\n", "xr = 1.0 # x-coordinate of the right boundary\n", "ncells = 100 # number of cells in the discretization\n", "nsteps = 100 # number of steps in the reactive transport simulation\n", "dx = (xr - xl) / ncells # length of the mesh cells (in units of m)\n", "dt = 30 * minute # time step\n", "\n", "# Physical parameters\n", "D = 1.0e-9 # diffusion coefficient (in units of m2/s)\n", "v = 1.0 / week # fluid pore velocity (in units of m/s)\n", "T = 60.0 + 273.15 # temperature (in units of K)\n", "P = 100 * 1e5 # pressure (in units of Pa)\n", "phi = 0.1 # the porosity" ] }, { "cell_type": "markdown", "id": "2f970197", "metadata": {}, "source": [ "Next, we generate the coordinates of the mesh nodes (array `xcells`) by equally dividing the interval *[xr, xl]* with\n", "the number of cells `ncells`. The length between each consecutive mesh node is computed and stored in `dx` (the\n", "length of the mesh cells)." ] }, { "cell_type": "code", "execution_count": null, "id": "e91d433b", "metadata": {}, "outputs": [], "source": [ "xcells = np.linspace(xl, xr, ncells + 1) # interval [xl, xr] split into ncells" ] }, { "cell_type": "markdown", "id": "b4a06199", "metadata": {}, "source": [ "The boolean variable `dirichlet` is set to `True` or `False` depending on which boundary condition is considered in\n", "the numerical calculation. `False` corresponds to imposing the flux of the injected fluid, otherwise, `True` means\n", "imposing the composition of the fluid on the left boundary." ] }, { "cell_type": "code", "execution_count": null, "id": "c451d572", "metadata": {}, "outputs": [], "source": [ "dirichlet = False # parameter that determines whether Dirichlet BC must be used" ] }, { "cell_type": "markdown", "id": "f434547c", "metadata": {}, "source": [ "To make sure that the applied finite-volume scheme is stable, we need to keep track of Courant–Friedrichs–Lewy (CFL)\n", "number, which should be less than 1.0." ] }, { "cell_type": "code", "execution_count": null, "id": "7ca27b06", "metadata": {}, "outputs": [], "source": [ "CFL = v * dt / dx\n", "assert CFL <= 1.0, f\"Make sure that CFL = {CFL} is less that 1.0\"" ] }, { "cell_type": "markdown", "id": "ad2a279f", "metadata": {}, "source": [ "## Specifying the quantities and properties to be outputted\n", "\n", "Before running the reactive transport simulations, we specify the list of parameters we are interested in\n", "outputting. In this case, it is pH, molality of `H+`, `Ca++`, `Mg++`, `HCO3-`, `CO2(aq)`, as well as a phase volume\n", "of calcite and dolomite." ] }, { "cell_type": "code", "execution_count": null, "id": "575c6ef8", "metadata": {}, "outputs": [], "source": [ "output_quantities = \"\"\"\n", " pH\n", " speciesMolality(H+)\n", " speciesMolality(Ca++)\n", " speciesMolality(Mg++)\n", " speciesMolality(HCO3-)\n", " speciesMolality(CO2(aq))\n", " phaseVolume(Calcite)\n", " phaseVolume(Dolomite)\n", " speciesMolality(CO3--)\n", " speciesMolality(CaCl+)\n", " speciesMolality(Ca(HCO3)+)\n", " speciesMolality(MgCl+)\n", " speciesMolality(Mg(HCO3)+)\n", " speciesMolality(OH-)\n", "\"\"\".split()" ] }, { "cell_type": "markdown", "id": "8fd34cc9", "metadata": {}, "source": [ "Then, we define the list of name for the DataFrame columns. Note, that they must correspond\n", "to the order of the properties defined in the `output_quantities` list:" ] }, { "cell_type": "code", "execution_count": null, "id": "ce05fd79", "metadata": {}, "outputs": [], "source": [ "column_quantities = \"\"\"\n", " pH\n", " Hcation\n", " Cacation\n", " Mgcation\n", " HCO3anion\n", " CO2aq\n", " calcite\n", " dolomite\n", " CO3anion\n", " CaClcation\n", " CaHCO3cation\n", " MgClcation\n", " MgHCO3cation\n", " OHanion\n", "\"\"\".split()" ] }, { "cell_type": "code", "execution_count": null, "id": "840c0582", "metadata": {}, "outputs": [], "source": [ "# Create the list of columns stored in dataframes\n", "columns = ['step', 'x'] + column_quantities\n", "import pandas as pd" ] }, { "cell_type": "code", "execution_count": null, "id": "641bbb75", "metadata": {}, "outputs": [], "source": [ "# Initialize dataframes with above defined columns\n", "df = pd.DataFrame(columns=columns)" ] }, { "cell_type": "markdown", "id": "40a7e4dd", "metadata": {}, "source": [ "## Organization of the program\n", "\n", "The main part of the program (at the bottom of this tutorial) consists of three parts, each represented by a Python\n", "function and documented in the following sections:\n", "* creation of folders for the results (function `make_results_folders()`),\n", "* simulation of reactive transport problem (method `simulate()`), and\n", "* plotting of the obtained results.\n", "\n", "## Creating folders for the outputted results\n", "\n", "Using **os** package, we create required folders for outputting the obtained results and for the plot and video\n", "files later." ] }, { "cell_type": "code", "execution_count": null, "id": "d8b8948b", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "folder_results = 'results-rt-calcite-dolomite'\n", "def make_results_folders():\n", " os.system('mkdir -p ' + folder_results)" ] }, { "cell_type": "markdown", "id": "d5f255f8", "metadata": {}, "source": [ "## Performing the reactive transport simulation\n", "\n", "The reactive transport simulation is performed in the function `simulate`, which consists from the several building\n", "blocks (functions):\n", "* initialization of the reactive transport problem and\n", "* performing the reactive transport simulation along defined time interval.\n", "\n", "The preparatory initialization step consists of the following sub-steps:\n", "* definition of chemical system with its phases and species using `define_chemical_system()`,\n", "* definition of the initial condition of the reactive transport problem in `define_initial_condition()`,\n", "* definition of the boundary condition of the reactive transport problem in `define_initial_condition()`,\n", "* generation of auxiliary indices to partition elements using `partition_indices()` and elements' partitioning \n", "corresponding to fluid and solid species with function `partition_elements_in_mesh_cell()`, and finally\n", "* definition of instance of [EquilibriumSolver](https://reaktoro.org/cpp/classReaktoro_1_1EquilibriumSolver.html).\n", "\n", "The simulation of the reactive transport problem is represented by the loop over discretized time interval until\n", "the final time is reached. On each step of this loop, the following functionality of performed:\n", "* transport calculations using method `transport()`,\n", "* reactive chemical calculations with `reactive_chemistry()` function, and\n", "* saving the obtained results by means of `outputstate()`." ] }, { "cell_type": "markdown", "id": "c7ff588a", "metadata": {}, "source": [ "Performing of the transport and reactive chemistry sequentially is possible due to the *operator splitting\n", "procedure*, in which we first update the amounts of elements `b`. These updated amounts of\n", "elements in the cell are used to evaluate its new chemical equilibrium state, thus producing new amounts of the\n", "species in both the fluid and solid phases (available in the list `states` of\n", "[ChemicalState](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html) objects). This chemical reaction\n", "equilibrium calculation step, at each mesh cell, permits, for example, aqueous species and minerals to react,\n", "and thus causes mineral dissolution or precipitation, depending on how much the amount of mineral species changes.\n", "This can then be used, for example, to compute new porosity value for the cell." ] }, { "cell_type": "code", "execution_count": null, "id": "b22f6345", "metadata": {}, "outputs": [], "source": [ "def simulate():\n", " # Construct the chemical system with its phases and species\n", " system = define_chemical_system()\n", "\n", " # Define the initial condition of the reactive transport modeling problem\n", " state_ic = define_initial_condition(system)\n", "\n", " # Define the boundary condition of the reactive transport modeling problem\n", " state_bc = define_boundary_condition(system)\n", "\n", " # Generate indices of partitioning fluid and solid species\n", " nelems, ifluid_species, isolid_species = partition_indices(system)\n", "\n", " # Partitioning fluid and solid species\n", " b, bfluid, bsolid, b_bc = partition_elements_in_mesh_cell(ncells, nelems, state_ic, state_bc)\n", "\n", " # Create a list of chemical states for the mesh cells (one for each cell, initialized to state_ic)\n", " states = [state_ic.clone() for _ in range(ncells + 1)]\n", "\n", " # Create the equilibrium solver object for the repeated equilibrium calculation\n", " solver = EquilibriumSolver(system)\n", "\n", " # Running the reactive transport simulation loop\n", " step = 0 # the current step number\n", " t = 0.0 # the current time (in seconds)\n", "\n", " # Output the initial state of the reactive transport calculation\n", " outputstate_df(step, system, states)\n", "\n", " with tqdm(total=nsteps, desc=\"Reactive transport simulations\") as pbar:\n", " while step <= nsteps:\n", " # Perform transport calculations\n", " bfluid, bsolid, b = transport(states, bfluid, bsolid, b, b_bc, nelems, ifluid_species, isolid_species)\n", "\n", " # Perform reactive chemical calculations\n", " states = reactive_chemistry(solver, states, b)\n", "\n", " # Increment time step and number of time steps\n", " t += dt\n", " step += 1\n", "\n", " # Output the current state of the reactive transport calculation\n", " outputstate_df(step, system, states)\n", "\n", " # Update a progress bar\n", " pbar.update(1)" ] }, { "cell_type": "markdown", "id": "e6093bc7", "metadata": {}, "source": [ "Subsections below correspond to the methods responsible for each of the functional parts of `simulate()` method.\n", "\n", "### Construction of the chemical system with its phases and species\n", "\n", "Reaktoro is a general-purpose chemical solver that avoids as much as possible presuming specific assumptions about\n", "your problems. Thus, you need to specify how your chemical system should be defined. This encompasses the\n", "specification of all phases in the system as well as the chemical species that compose each phase. By using the\n", "[ChemicalEditor](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalEditor.html) class, you can conveniently achieve\n", "this as shown below in method `define_chemical_system()`.\n", "\n", "In this step, we create an object of class [ChemicalEditor](\n", "https://reaktoro.org/cpp/classReaktoro_1_1ChemicalEditor.html) and specify two phases, an *aqueous* and a\n", "*mineral*, should be considered in the chemical system. The aqueous phase is defined by using a list of\n", "compounds, which will be broken into a list of element names, and the database will then be searched for all\n", "species that could be formed out of those elements. The mineral phase is defined as three mineral species:\n", "quartz (SiO2), calcite (CaCO3), and dolomite (CaMg(CO3)2).\n", "\n", "> **Note**: An automatic search for chemical species can result in a large number of species in the phase,\n", "> potentially causing the chemical reaction calculations to be more computationally expensive. If you are using\n", "> Reaktoro in demanding applications (e.g., as a chemical solver in a reactive transport simulator), you might want\n", "> to manually specify the chemical species of each phase in your chemical system. This can be achieved by providing a\n", "> list of species names as `editor.addAqueousPhaseWithElements(['H2O(l)', 'H+', 'OH-', 'Na+', 'Cl-', 'Ca++', 'Mg++',\n", "> 'HCO3-', 'CO2(aq)', 'CO3--' ]).\n", "\n", "Finally, we create an object of class [ChemicalSystem](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html)\n", "using the chemical system definition details stored in the object `editor`.\n", "\n", "> **Note**: [ChemicalSystem](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html) is perhaps the main class\n", "> in Reaktoro. An object of this class stores the phases, species, and elements in our defined chemical system,\n", "> as well as provides the means to compute many types of thermodynamic properties, such as *standard thermodynamic\n", "> properties* (e.g., standard Gibbs energies, standard enthalpies, and standard volumes of species),\n", "> and *thermo-chemical properties* (e.g., activity and activity coefficients of species; density, enthalpy and\n", "> internal energy of phases). As you learn more about other Reaktoro's classes, you will note that an object of class\n", "> [ChemicalSystem](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html) is almost always needed for their\n", "> initialization.\n", "\n", "The *activity coefficients* of the aqueous species in this tutorial are calculated using the\n", "*Pitzer model* (unlike the default *HKF extended Debye-Huckel model*) for solvent water and ionic species,\n", "except for the aqueous species CO2(aq), for which the *Drummond model* is used.\n", "The *standard chemical potentials* of the species are calculated using the equations of state of Helgeson and\n", "Kirkham (1974), Helgeson et al. (1974), Tanger and Helgeson (1988), Shock and Helgeson (1988), and Shock et al. (\n", "1992). The database file [slop98.dat](https://github.com/reaktoro/reaktoro/blob/master/databases/supcrt/slop98.dat)\n", "from the software SUPCRT92 is used to obtain the parameters for the equations of state.\n", "The equation of state of Wagner and Pruss (2002) is used to calculate the *density of water* and its temperature and\n", "pressure derivatives. Kinetics of *dissolution* and *precipitation* of both calcite and dolomite is neglected, i.e.,\n", "the local equilibrium assumption is employed." ] }, { "cell_type": "code", "execution_count": null, "id": "6b3806b5", "metadata": {}, "outputs": [], "source": [ "def define_chemical_system():\n", " # Construct the chemical system with its phases and species\n", " db = Database('supcrt98.xml')\n", "\n", " editor = ChemicalEditor(db)\n", "\n", " editor.addAqueousPhaseWithElements('H O Na Cl Ca Mg C Si Ca') \\\n", " .setChemicalModelPitzerHMW() \\\n", " .setActivityModelDrummondCO2()\n", " editor.addMineralPhase('Quartz')\n", " editor.addMineralPhase('Calcite')\n", " editor.addMineralPhase('Dolomite')\n", "\n", " system = ChemicalSystem(editor)\n", "\n", " return system" ] }, { "cell_type": "markdown", "id": "c9cbcde5", "metadata": {}, "source": [ "### Initial condition (IC) of the reactive transport problem\n", "\n", "We have now defined and constructed our chemical system of interest, enabling us to move on to the next step in\n", "Reaktoro's modeling workflow: *defining our chemical reaction problems*. Below we define its **initial condition**\n", "with already prescribed equilibrium conditions for *temperature*, *pressure*, and *amounts of elements* that are\n", "consistent to model reactive transport of injected NaCl-MgCl2-CaCl2 brine into the\n", "rock-fluid composition of quartz and calcite at 60 °C and 100 bar. In particular, we consider resident fluid is a\n", "0.7 molal NaCl brine in equilibrium with the rock minerals with a calculated pH of 8.43095." ] }, { "cell_type": "code", "execution_count": null, "id": "445c8178", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def define_initial_condition(system):\n", " problem_ic = EquilibriumProblem(system)\n", " problem_ic.setTemperature(T)\n", " problem_ic.setPressure(P)\n", " problem_ic.add('H2O', 1.0, 'kg')\n", " problem_ic.add('NaCl', 0.7, 'mol')\n", " problem_ic.add('CaCO3', 10, 'mol')\n", " problem_ic.add('SiO2', 10, 'mol')\n", "\n", " # Calculate the equilibrium states for the initial conditions\n", " state_ic = equilibrate(problem_ic)\n", "\n", " # Scale the volumes of the phases in the initial condition\n", " state_ic.scalePhaseVolume('Aqueous', 0.1, 'm3')\n", " state_ic.scalePhaseVolume('Quartz', 0.882, 'm3')\n", " state_ic.scalePhaseVolume('Calcite', 0.018, 'm3')\n", "\n", " return state_ic" ] }, { "cell_type": "markdown", "id": "4749a886", "metadata": {}, "source": [ "> **Note**: Our **0.7 molal NaCl brine** is here represented by the mixture of 1 kg of H2O and 0.7 mol of\n", "> NaCl.\n", "\n", "> **Note**: After providing the amounts of substances H2O, NaCl, quartz SiO2,\n", "> and calcite CaCO3 in the above code, Reaktoro parses these chemical formulas and determines the\n", "> elements and their coefficients. Once this is done, the amount of each element stored inside the object\n", "> `problem_ic` is incremented according to the given amount of substance and its coefficient in the formula. The\n", "> amounts of elements you provide are then used as constraints for the Gibbs energy minimization calculation when\n", "> computing the state of chemical equilibrium (i.e., when we try to find the amounts of all species in the system\n", "> that corresponds to a state of minimum Gibbs energy and at the same time satisfying the *element amounts\n", "> constraints*).\n", "\n", "> **Note**: Please note that we are not condemning the input form shown above in terms of element amounts,\n", "> but only telling you to be attentive with the values you input. If you are using Reaktoro as a chemical reaction\n", "> solver in a reactive transport simulator, for example, you'll most likely need to work directly with given amounts\n", "> of elements, which shows that this input form is required in certain cases. For such time-dependent modeling\n", "> problems, you often only need to ensure that the initial conditions for elements amounts result in feasible initial\n", "> species amounts.\n", "\n", "Next, we use method [equilibrate](https://reaktoro.org/cpp/namespaceReaktoro.html#af2d3b39d3e0b8f9cb5a4d9bbb06b697e)\n", "to calculate the chemical equilibrium state of the system with the given initial conditions stored in the objects\n", "`problem_ic`. The numerical solution of each problem results in the objects `state_ic` of class\n", "[ChemicalState](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html), which stores the temperature,\n", "pressure, and the amounts of every species in the system.\n", "For this calculation, Reaktoro uses an efficient **Gibbs energy minimization** computation to determine the species\n", "amounts that correspond to a state of minimum Gibbs energy in the system, while satisfying the prescribed amount\n", "conditions for temperature, pressure, and element amounts.\n", "\n", "The function ends with scaling the volumes of the aqueous and mineral phases so that they are consistent with a 10 %\n", "porosity and the required volume percentages of the rock minerals (98 %vol of quartz and 2 %vol\n", "of calcite).\n", "\n", "### Boundary condition (BC) of the reactive transport problem\n", "\n", "Next, we define the **boundary condition** of the constructed chemical system with its *temperature*, *pressure*,\n", "and *amounts of elements*. In particular, we prescribe the amount of injected aqueous fluid resulting from\n", "the mixture of 1 kg of water with 0.90 moles of NaCl, 0.05 molal MgCl2, 0.01 CaCl2,\n", "and 0.75 molal CO2, in a state very close to CO2 saturation.\n", "The temperature and the pressure stays the same, i.e., 60 °C and 100 bar, respectively.\n", "\n", "After equilibration, the obtained chemical state representing the boundary condition for the injected fluid\n", "composition, we scale its volume to 1 m3. This is done so that the amounts of the species in the fluid are\n", "consistent with a \\mathrm{mol/m^3} scale." ] }, { "cell_type": "code", "execution_count": null, "id": "f6f8f07e", "metadata": {}, "outputs": [], "source": [ "def define_boundary_condition(system):\n", " # Define the boundary condition of the reactive transport modeling problem\n", " problem_bc = EquilibriumProblem(system)\n", " problem_bc.setTemperature(T)\n", " problem_bc.setPressure(P)\n", " problem_bc.add('H2O', 1.0, 'kg')\n", " problem_bc.add('NaCl', 0.90, 'mol')\n", " problem_bc.add('MgCl2', 0.05, 'mol')\n", " problem_bc.add('CaCl2', 0.01, 'mol')\n", " problem_bc.add('CO2', 0.75, 'mol')\n", "\n", " # Calculate the equilibrium states for the boundary conditions\n", " state_bc = equilibrate(problem_bc)\n", " # Scale the boundary condition state to 1 m3\n", " state_bc.scaleVolume(1.0, 'm3')\n", "\n", " return state_bc" ] }, { "cell_type": "markdown", "id": "aeca3075", "metadata": {}, "source": [ "### Indices of partitioning fluid and solid species\n", "\n", "Only species in fluid phases are mobile and transported by advection and diffusion mechanisms. The solid phases are\n", "immobile. The code below identifies the indices of the fluid and solid species. We use methods\n", "[indicesFluidSpecies](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html#ac2a8b713f46f7a66b2731ba63faa95ad)\n", "and [indicesSolidSpecies](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html#a8b0c237fff1d827f7bf2dbc911fa5bbf)\n", "of class [ChemicalSystem](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalSystem.html) to get the indices of the\n", "fluid and solid species, which are stored in the lists `ifluid_species` and `isolid_species`, respectively." ] }, { "cell_type": "code", "execution_count": null, "id": "d19a5856", "metadata": {}, "outputs": [], "source": [ "def partition_indices(system):\n", " nelems = system.numElements()\n", "\n", " ifluid_species = system.indicesFluidSpecies()\n", " isolid_species = system.indicesSolidSpecies()\n", "\n", " return nelems, ifluid_species, isolid_species" ] }, { "cell_type": "markdown", "id": "35e14c57", "metadata": {}, "source": [ "### Partitioning fluid and solid species\n", "\n", "In this function, we create arrays to keep track of the amounts of elements in the fluid and solid partition\n", "(i.e., the amounts of elements among all fluid phases, here only an aqueous phase, and the amounts of elements among\n", "all solid phases, here the mineral phases). For that, we define the arrays `b`, `bfluid`, `bsolid`, that\n", "will store, respectively, the concentrations (mol/m3) of each element in the system, in the fluid\n", "partition, and in the solid partition at every time step.\n", "\n", "The array `b` is initialized with the concentrations of the elements at the initial chemical state, `state_ic`,\n", "using method [elementAmounts](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html#a827457e68a90f89920c13f0cc06fda78)\n", "of class [ChemicalState](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html). The array `b_bc` stores\n", "the concentrations of each element on the boundary in mol/m3fluid." ] }, { "cell_type": "code", "execution_count": null, "id": "36002a59", "metadata": {}, "outputs": [], "source": [ "def partition_elements_in_mesh_cell(ncells, nelems, state_ic, state_bc):\n", " # The concentrations of each element in each mesh cell (in the current time step)\n", " b = np.zeros((ncells, nelems))\n", " # Initialize the concentrations (mol/m3) of the elements in each mesh cell\n", " b[:] = state_ic.elementAmounts()\n", "\n", " # The concentrations (mol/m3) of each element in the fluid partition, in each mesh cell\n", " bfluid = np.zeros((ncells, nelems))\n", "\n", " # The concentrations (mol/m3) of each element in the solid partition, in each mesh cell\n", " bsolid = np.zeros((ncells, nelems))\n", "\n", " # Initialize the concentrations (mol/m3) of each element on the boundary\n", " b_bc = state_bc.elementAmounts()\n", "\n", " return b, bfluid, bsolid, b_bc" ] }, { "cell_type": "markdown", "id": "f20304a3", "metadata": {}, "source": [ "### Reactive transport cycle\n", "\n", "#### Transport\n", "\n", "This step updates in the fluid partition `bfluid` using the transport equations (without reactions).\n", "The `transport_fullimplicit()` function below is responsible for solving an advection-diffusion equation, that is\n", "later applied to transport the concentrations (mol/m3) of elements in the fluid partition (*a\n", "simplification that is possible because of common diffusion coefficients and velocities of the fluid species,\n", "otherwise the transport of individual fluid species would be needed*).\n", "\n", "To match the units of concentrations of the elements in the fluid measure in mol/m3bulk and the\n", "imposed concentration `b_bc[j]` mol/m3fluid, e need to multiply it by the porosity `phi_bc` on the\n", "boundary cell m3fluid/m3bulk. We use function\n", "[properties](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html#ad3fa8fd9e1b948da7a698eb020513f3d)\n", "of the class [ChemicalState](https://reaktoro.org/cpp/classReaktoro_1_1ChemicalState.html) to retrieve fluid volume\n", "m3fluid and total volume m3bulk in the inflow boundary cell.\n", "\n", "The updated amounts of elements in the fluid partition are then summed with the amounts of elements in the solid partition\n", "`bsolid`, which remained constant during the transport step), and thus updating the amounts of elements in the\n", "chemical system `b`. Reactive transport calculations involve the solution of a system of\n", "advection-diffusion-reaction equations." ] }, { "cell_type": "code", "execution_count": null, "id": "952113b1", "metadata": {}, "outputs": [], "source": [ "def transport(states, bfluid, bsolid, b, b_bc, nelems, ifluid_species, isolid_species):\n", " # Collect the amounts of elements from fluid and solid partitions\n", " for icell in range(ncells):\n", " bfluid[icell] = states[icell].elementAmountsInSpecies(ifluid_species)\n", " bsolid[icell] = states[icell].elementAmountsInSpecies(isolid_species)\n", "\n", " # Get the porosity of the boundary cell\n", " bc_cell = 0\n", " phi_bc = states[bc_cell].properties().fluidVolume().val / states[bc_cell].properties().volume().val\n", "\n", " # Transport each element in the fluid phase\n", " for j in range(nelems):\n", " transport_fullimplicit(bfluid[:, j], dt, dx, v, D, phi_bc * b_bc[j])\n", "\n", " # Update the amounts of elements in both fluid and solid partitions\n", " b[:] = bsolid + bfluid\n", "\n", " return bfluid, bsolid, b" ] }, { "cell_type": "markdown", "id": "aaabfdb6", "metadata": {}, "source": [ "##### Transport calculation with finite-volume scheme\n", "\n", "The function `transport()` expects a conservative property (argument `u`) (e.g., the concentration mol/m3\n", "of *j*th element in the fluid given by `bfluid[j]`), the time step (`dt`), the mesh cell length (`dx`),\n", "the fluid velocity (`v`), the diffusion coefficient (`D`), and the boundary condition of the conservative property\n", "(`g`) (e.g., the concentration of the *j*th element in the fluid on the left boundary).\n", "\n", "The transport equations are solved with a finite volume method, where diffusion and convection are treated implicitly.\n", "Its discretization in space and time (implicit) results in the constants `alpha` and `beta`. These correspond to\n", "the diffusion and advection terms in the equation: `D*dt/dx**2` and `v*dt/dx`, respectively.\n", "\n", "Arrays `a`, `b`, `c` are the diagonals in the tridiagonal matrix that results by writing all discretized equations\n", "in a matrix equation. This system of linear equations is solved by the tridiagonal matrix algorithm, also known\n", "as the Thomas algorithm." ] }, { "cell_type": "code", "execution_count": null, "id": "14de3192", "metadata": {}, "outputs": [], "source": [ "def transport_fullimplicit(u, dt, dx, v, D, ul):\n", " # Number of DOFs\n", " n = len(u)\n", " alpha = D * dt / dx ** 2\n", " beta = v * dt / dx\n", "\n", " # Upwind finite volume scheme\n", " a = np.full(n, -beta - alpha)\n", " b = np.full(n, 1 + beta + 2 * alpha)\n", " c = np.full(n, -alpha)\n", "\n", " # Set the boundary condition on the left cell\n", " if dirichlet:\n", " # Use Dirichlet BC boundary conditions\n", " b[0] = 1.0\n", " c[0] = 0.0\n", " u[0] = ul\n", "\n", " else:\n", " # Flux boundary conditions (implicit scheme for the advection)\n", " # Left boundary\n", " b[0] = 1 + alpha + beta\n", " c[0] = -alpha # stays the same as it is defined -alpha\n", " u[0] += beta * ul # = dt/dx * v * g, flux that we prescribe is equal v * ul\n", "\n", " # Right boundary is free\n", " a[-1] = - beta\n", " b[-1] = 1 + beta\n", "\n", " # Solve a tridiagonal matrix equation\n", " thomas(a, b, c, u)" ] }, { "cell_type": "markdown", "id": "f9095759", "metadata": {}, "source": [ "##### Solving the system of equations obtained from finite volume discretization\n", "\n", "The tridiagonal matrix equation is solved using the Thomas algorithm (or the TriDiagonal Matrix Algorithm (TDMA)).\n", "It is a simplified form of Gaussian elimination that can be used to solve tridiagonal systems of equations." ] }, { "cell_type": "code", "execution_count": null, "id": "425da0d4", "metadata": {}, "outputs": [], "source": [ "def thomas(a, b, c, d):\n", " n = len(d)\n", " c[0] /= b[0]\n", " for i in range(1, n - 1):\n", " c[i] /= b[i] - a[i] * c[i - 1]\n", " d[0] /= b[0]\n", " for i in range(1, n):\n", " d[i] = (d[i] - a[i] * d[i - 1]) / (b[i] - a[i] * c[i - 1])\n", " x = d\n", " for i in reversed(range(0, n - 1)):\n", " x[i] -= c[i] * x[i + 1]\n", " return x" ] }, { "cell_type": "markdown", "id": "0a33ddd0", "metadata": {}, "source": [ "#### Reactive chemistry\n", "\n", "The chemical equilibrium calculations performed in each mesh cell, using *Gibbs energy minimisation* algorithm (\n", "provided by the class [EquilibriumSolver](https://reaktoro.org/cpp/classReaktoro_1_1EquilibriumSolver.html))." ] }, { "cell_type": "code", "execution_count": null, "id": "da74b74b", "metadata": {}, "outputs": [], "source": [ "def reactive_chemistry(solver, states, b):\n", " # Equilibrating all cells with the updated element amounts\n", " for icell in range(ncells):\n", " solver.solve(states[icell], T, P, b[icell])\n", " return states" ] }, { "cell_type": "markdown", "id": "f14de43a", "metadata": {}, "source": [ "### Results saving and analyzing\n", "\n", "Function `outputstate_df` is the auxiliary function to add data to the DataFrame at each time step." ] }, { "cell_type": "code", "execution_count": null, "id": "ce74d971", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def outputstate_df(step, system, states):\n", " # Define the instance of ChemicalQuantity class\n", " quantity = ChemicalQuantity(system)\n", "\n", " # Create the list with empty values to populate with chemical properties\n", " values = [None] * len(columns)\n", " for state, x in zip(states, xcells):\n", "\n", " # Populate values with number of reactive transport step and spacial coordinates\n", " values[0] = step\n", " values[1] = x\n", "\n", " # Update the\n", " quantity.update(state)\n", " for quantity_name, i in zip(output_quantities, range(2, len(states))):\n", " values[i] = quantity.value(quantity_name) * (100 / (1 - phi) if \"phaseVolume\" in quantity_name else 1)\n", " df.loc[len(df)] = values" ] }, { "cell_type": "markdown", "id": "558de1bb", "metadata": {}, "source": [ "### Plotting of the obtained results\n", "\n", "The last block of the main routine is dedicated to the plotting of the results in a Jupyter app generated by the\n", "library **bokeh**. It is an interactive visualization library that provides elegant, concise construction of\n", "versatile graphics, and affords high-performance interactivity over large or streaming datasets.\n", "\n", "Below, we list auxiliary functions that we use in plotting. Function `titlestr` returns a string for the title\n", "of a figure in the format Time: #h##m" ] }, { "cell_type": "code", "execution_count": null, "id": "111744e9", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def titlestr(t):\n", " t = t / minute # Convert from seconds to minutes\n", " h = int(t) / 60 # The number of hours\n", " m = int(t) % 60 # The number of remaining minutes\n", " return 'Time: %2dh %2dm' % (h, m)" ] }, { "cell_type": "markdown", "id": "4544ecaa", "metadata": {}, "source": [ "Routines `plot_figures_ph()`, `plot_figures_calcite_dolomite()`, and 'plot_figures_aqueous_species()'\n", "are dedicated to drawing the plots with chemical properties on the selected steps that are specified by the user\n", "below." ] }, { "cell_type": "code", "execution_count": null, "id": "367259d1", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def plot_figures_ph(steps, files):\n", " # Plot ph on the selected steps\n", " plots = []\n", " for i in steps:\n", " t = i * dt\n", " source = ColumnDataSource(df[df['step'] == i])\n", "\n", " p = figure(plot_width=600, plot_height=250)\n", " p.line(source.data['x'], source.data['pH'], color='teal', line_width=2, legend_label='pH')\n", " p.x_range = Range1d(-0.001, 1.001)\n", " p.y_range = Range1d(2.5, 12.0)\n", " p.xaxis.axis_label = 'Distance [m]'\n", " p.yaxis.axis_label = 'pH'\n", " p.legend.location = 'bottom_right'\n", " p.title.text = titlestr(t)\n", "\n", " plots.append([p])\n", "\n", " grid = gridplot(plots)\n", " show(grid)" ] }, { "cell_type": "code", "execution_count": null, "id": "758ca94f", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def plot_figures_calcite_dolomite(steps, files):\n", " plots = []\n", " for i in steps:\n", " t = i * dt\n", " source = ColumnDataSource(df[df['step'] == i])\n", "\n", " p = figure(plot_width=600, plot_height=250)\n", " p.line(source.data['x'], source.data['calcite'], color='blue', line_width=2, legend_label='Calcite',\n", " muted_color='blue', muted_alpha=0.2)\n", " p.line(source.data['x'], source.data['dolomite'], color='orange', line_width=2, legend_label='Dolomite',\n", " muted_color='orange', muted_alpha=0.2)\n", " p.x_range = Range1d(-0.001, 1.001)\n", " p.y_range = Range1d(-0.1, 2.1)\n", " p.xaxis.axis_label = 'Distance [m]'\n", " p.yaxis.axis_label = 'Mineral Volume [%vol]'\n", " p.legend.location = 'center_right'\n", " p.title.text = titlestr(t)\n", " p.legend.click_policy = 'mute'\n", " plots.append([p])\n", "\n", " grid = gridplot(plots)\n", " show(grid)" ] }, { "cell_type": "code", "execution_count": null, "id": "f7ce9cd2", "metadata": {}, "outputs": [], "source": [ "def plot_figures_aqueous_species(steps, files):\n", " plots = []\n", " for i in steps:\n", " t = i * dt\n", " source = ColumnDataSource(df[df['step'] == i])\n", " p = figure(plot_width=600, plot_height=300, y_axis_type = 'log',)\n", "\n", " p.line(source.data['x'], source.data['Cacation'], color='blue', line_width=2, legend_label='Ca++')\n", " p.line(source.data['x'], source.data['Mgcation'], color='orange', line_width=2, legend_label='Mg++')\n", " p.line(source.data['x'], source.data['HCO3anion'], color='green', line_width=2, legend_label='HCO3-')\n", " p.line(source.data['x'], source.data['CO2aq'], color='red', line_width=2, legend_label='CO2(aq)')\n", " p.line(source.data['x'], source.data['Hcation'], color='darkviolet', line_width=2, legend_label='H+')\n", " p.line(source.data['x'], source.data['CO3anion'], color='teal', line_width=2, legend_label='CO3--')\n", " p.line(source.data['x'], source.data['CaClcation'], color='coral', line_width=2, legend_label='CaCl+')\n", " p.line(source.data['x'], source.data['CaHCO3cation'], color='mediumseagreen', line_width=2, legend_label='Ca(HCO3)+')\n", " p.line(source.data['x'], source.data['MgClcation'], color='darkred', line_width=2, legend_label='MgCl+')\n", " p.line(source.data['x'], source.data['MgHCO3cation'], color='indigo', line_width=2, legend_label='Mg(HCO3)+')\n", " p.line(source.data['x'], source.data['OHanion'], color='grey', line_width=2, legend_label='OH-')\n", "\n", " p.x_range = Range1d(-0.001, 1.001)\n", " p.y_range = Range1d(1e-9, 1e0)\n", " p.xaxis.axis_label = 'Distance [m]'\n", " p.yaxis.axis_label = 'Concentration [molal]'\n", " p.legend.location = 'top_right'\n", " p.title.text = titlestr(t)\n", " p.legend.click_policy = 'mute'\n", " plots.append([p])\n", "\n", " grid = gridplot(plots)\n", " show(grid)" ] }, { "cell_type": "markdown", "id": "26e6b904", "metadata": {}, "source": [ "# Main parts of the tutorial\n", "\n", "First, we create folders for the results:" ] }, { "cell_type": "code", "execution_count": null, "id": "818cf2e4", "metadata": {}, "outputs": [], "source": [ "make_results_folders()" ] }, { "cell_type": "markdown", "id": "e29243cc", "metadata": {}, "source": [ "Run the reactive transport simulations. Note that due to selected Pitzer activity model, the reactive transport\n", "simulation might be slower than, for instance, Debye-Huckel one." ] }, { "cell_type": "code", "execution_count": null, "id": "45e59947", "metadata": {}, "outputs": [], "source": [ "simulate()" ] }, { "cell_type": "markdown", "id": "91f2d9c9", "metadata": {}, "source": [ "To inspect the collected values on any reactive transport step, one can run the code below, where we first define\n", "the number of the step and then select the columns of DataFrame we are interested in (cells indices together with the\n", "chemical properties):" ] }, { "cell_type": "code", "execution_count": null, "id": "7cfa1df1", "metadata": {}, "outputs": [], "source": [ "step = 0\n", "df_step = df[df['step'] == step].loc[:, ['x'] + column_quantities]\n", "df_step" ] }, { "cell_type": "markdown", "id": "9d3658db", "metadata": {}, "source": [ "Select the steps, on which results must plotted:" ] }, { "cell_type": "code", "execution_count": null, "id": "088e7924", "metadata": {}, "outputs": [], "source": [ "selected_steps_to_plot = [10, 100]\n", "assert all(step <= nsteps for step in selected_steps_to_plot), f\"Make sure that selceted steps are less than \" \\\n", " f\"total amount of steps {nsteps}\"" ] }, { "cell_type": "markdown", "id": "0a3b5083", "metadata": {}, "source": [ "Then, we collect files with results corresponding to each step:" ] }, { "cell_type": "code", "execution_count": null, "id": "602eca75", "metadata": {}, "outputs": [], "source": [ "print(\"Collecting files...\")\n", "files = [file for file in natsorted(os.listdir(folder_results))]" ] }, { "cell_type": "markdown", "id": "15d2ffcf", "metadata": {}, "source": [ "Outputting the plots to the notebook requires the call of `output_notebook()` that specifies outputting the plot\n", "inline in the Jupyter notebook:" ] }, { "cell_type": "code", "execution_count": null, "id": "453bc86c", "metadata": {}, "outputs": [], "source": [ "output_notebook()" ] }, { "cell_type": "markdown", "id": "c2de5814", "metadata": {}, "source": [ "Plot ph on the selected steps:" ] }, { "cell_type": "code", "execution_count": null, "id": "ef49b933", "metadata": {}, "outputs": [], "source": [ "plot_figures_ph(selected_steps_to_plot, files)" ] }, { "cell_type": "markdown", "id": "64cbf5cb", "metadata": {}, "source": [ "Plot calcite and dolomite on the selected steps:" ] }, { "cell_type": "code", "execution_count": null, "id": "597cdd71", "metadata": {}, "outputs": [], "source": [ "plot_figures_calcite_dolomite(selected_steps_to_plot, files)" ] }, { "cell_type": "markdown", "id": "c961bfa6", "metadata": {}, "source": [ "Plot aqueous species on the selected steps:" ] }, { "cell_type": "code", "execution_count": null, "id": "a2ae4946", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "plot_figures_aqueous_species(selected_steps_to_plot, files)" ] }, { "cell_type": "markdown", "id": "f80923fc", "metadata": {}, "source": [ "To study the time-dependent behavior of the chemical properties, we create a Bokeh application using function\n", "`modify_doc(doc)`. It creates Bokeh content and adds it to the app. Streaming of the reactive transport data:" ] }, { "cell_type": "code", "execution_count": null, "id": "52050dc4", "metadata": { "lines_to_next_cell": 1 }, "outputs": [], "source": [ "def modify_doc(doc):\n", " # Initialize the data by the initial chemical state\n", " # source = ColumnDataSource(df[df['step'] == 0].loc[:, column_quantities])\n", " source = ColumnDataSource(df[df['step'] == 0])\n", "\n", " # Auxiliary function that returns a string for the title of a figure in the format Time: #h##m\n", " def titlestr(t):\n", " t = t / minute # Convert from seconds to minutes\n", " h = int(t) / 60 # The number of hours\n", " m = int(t) % 60 # The number of remaining minutes\n", " return 'Time: %2dh %2dm' % (h, m)\n", "\n", " # Plot for ph\n", " p1 = figure(plot_width=600, plot_height=250)\n", " p1.line(x='x', y='pH',\n", " color='teal', line_width=2, legend_label='pH',\n", " source=source)\n", " p1.x_range = Range1d(-0.001, 1.001)\n", " p1.y_range = Range1d(4.0, 9.0)\n", " p1.xaxis.axis_label = 'Distance [m]'\n", " p1.yaxis.axis_label = 'pH'\n", " p1.legend.location = 'bottom_right'\n", " p1.title.text = titlestr(0 * dt)\n", "\n", " # Plot for calcite and dolomite\n", " p2 = figure(plot_width=600, plot_height=250)\n", " p2.line(x='x', y='calcite', color='blue', line_width=2,\n", " legend_label='Calcite',\n", " muted_color='blue', muted_alpha=0.2,\n", " source=source)\n", " p2.line(x='x', y='dolomite', color='orange', line_width=2,\n", " legend_label='Dolomite',\n", " muted_color='orange', muted_alpha=0.2,\n", " source=source)\n", " p2.x_range = Range1d(-0.001, 1.001)\n", " p2.y_range = Range1d(-0.1, 2.1)\n", " p2.xaxis.axis_label = 'Distance [m]'\n", " p2.yaxis.axis_label = 'Mineral Volume [%vol]'\n", " p2.legend.location = 'center_right'\n", " p2.title.text = titlestr(0 * dt)\n", " p2.legend.click_policy = 'mute'\n", "\n", " # Plot for aqueous species\n", " p3 = figure(plot_width=600, plot_height=300, y_axis_type='log')\n", " p3.line(x='x', y='Cacation', color='blue', line_width=2, legend_label='Ca++', source=source)\n", " p3.line(x='x', y='Mgcation', color='orange', line_width=2, legend_label='Mg++', source=source)\n", " p3.line(x='x', y='HCO3anion', color='green', line_width=2, legend_label='HCO3-', source=source)\n", " p3.line(x='x', y='CO2aq', color='red', line_width=2, legend_label='CO2(aq)', source=source)\n", " p3.line(x='x', y='Hcation', color='darkviolet', line_width=2, legend_label='H+', source=source)\n", " p3.line(x='x', y='CO3anion', color='teal', line_width=2, legend_label='CO3--', source=source)\n", " p3.line(x='x', y='CaClcation', color='coral', line_width=2, legend_label='CaCl+', source=source)\n", " p3.line(x='x', y='CaHCO3cation', color='mediumseagreen', line_width=2, legend_label='Ca(HCO3)+', source=source)\n", " p3.line(x='x', y='MgClcation', color='darkred', line_width=2, legend_label='MgCl+', source=source)\n", " p3.line(x='x', y='MgHCO3cation', color='indigo', line_width=2, legend_label='Mg(HCO3)+', source=source)\n", " p3.line(x='x', y='OHanion', color='grey', line_width=2, legend_label='OH-', source=source)\n", "\n", " p3.x_range = Range1d(-0.001, 1.001)\n", " p3.y_range = Range1d(1e-9, 1e0)\n", " p3.xaxis.axis_label = 'Distance [m]'\n", " p3.yaxis.axis_label = 'Concentration [molal]'\n", " p3.legend.location = 'top_right'\n", " p3.title.text = titlestr(0 * dt)\n", " p3.legend.click_policy = 'mute'\n", "\n", " layout = column(p1, p2, p3)\n", "\n", " # Function that return the data dictionary with provided index of the file\n", " def update():\n", "\n", " if source.data['step'][0] + 1 <= nsteps:\n", " step_number = source.data['step'][0] + 1\n", " else:\n", " step_number = 0\n", "\n", " new_source = ColumnDataSource(df[df['step'] == step_number])\n", " new_data = dict(index=np.linspace(0, ncells, ncells + 1, dtype=int),\n", " step=new_source.data['step'],\n", " x=new_source.data['x'],\n", " pH=new_source.data['pH'],\n", " calcite=new_source.data['calcite'],\n", " dolomite=new_source.data['dolomite'],\n", " Hcation=new_source.data['Hcation'],\n", " Cacation=new_source.data['Cacation'],\n", " Mgcation=new_source.data['Mgcation'],\n", " HCO3anion=new_source.data['HCO3anion'],\n", " CO2aq=new_source.data['CO2aq'],\n", " CO3anion=new_source.data['CO3anion'],\n", " CaClcation=new_source.data['CaClcation'],\n", " CaHCO3cation=new_source.data['CaHCO3cation'],\n", " MgClcation=new_source.data['MgClcation'],\n", " MgHCO3cation=new_source.data['MgHCO3cation'],\n", " OHanion=new_source.data['OHanion'])\n", "\n", " p1.title.text = titlestr(step_number * dt)\n", " p2.title.text = titlestr(step_number * dt)\n", " p3.title.text = titlestr(step_number * dt)\n", "\n", " source.stream(new_data, rollover=ncells+1)\n", "\n", " doc.add_periodic_callback(update, 500)\n", " doc.add_root(layout)" ] }, { "cell_type": "markdown", "id": "46c10016", "metadata": {}, "source": [ "Outputting the plots to the notebook requires the call of `output_notebook()` that specifies outputting the plot\n", "inline in the Jupyter notebook. Finally, the function `modify_doc()` must be passed to `show`, so that the app defined\n", "by it is displayed inline.\n", "\n", "> **Important:** If you run this tutorial in the *localhost*, make sure that number provided to the variable\n", "`notebook_url` below coincides with the number of the localhost you have in your browser.\n", "\n", "In the app below, we refresh the reactive time step in a loop, which automatically updates the data source for the\n", "plots for ph, volume phases of calcite and dolomite, and mollalities of aqueous species (in logarithmic scale)." ] }, { "cell_type": "code", "execution_count": null, "id": "82cb5c50", "metadata": {}, "outputs": [], "source": [ "output_notebook()\n", "show(modify_doc, notebook_url=\"http://localhost:8888\")" ] } ], "metadata": { "kernelspec": { "display_name": "Python 3", "language": "python", "name": "python3" } }, "nbformat": 4, "nbformat_minor": 5 }