.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/howto_conversions/plot_B_feflowlib_BC_mesh.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_auto_examples_howto_conversions_plot_B_feflowlib_BC_mesh.py: Feflowlib: How to modify boundary conditions after conversion of a FEFLOW model. ================================================================================ .. sectionauthor:: Julian Heinze (Helmholtz Centre for Environmental Research GmbH - UFZ) This example shows how boundary conditions can be modified after converting a FEFLOW model. First we will change the values of the boundary conditions and later we will show how to define a new boundary mesh. .. GENERATED FROM PYTHON SOURCE LINES 12-13 0. Necessary imports .. GENERATED FROM PYTHON SOURCE LINES 13-23 .. code-block:: Python import tempfile from pathlib import Path import numpy as np import pyvista as pv import ogstools as ot from ogstools.examples import feflow_model_2D_HT from ogstools.feflowlib import assign_bulk_ids .. GENERATED FROM PYTHON SOURCE LINES 24-26 1. Load a FEFLOW model (.fem) as a FeflowModel object to further work it. During the initialisation, the FEFLOW file is converted. .. GENERATED FROM PYTHON SOURCE LINES 26-31 .. code-block:: Python temp_dir = Path(tempfile.mkdtemp("converted_models")) feflow_model = ot.FeflowModel(feflow_model_2D_HT, temp_dir / "HT_model") mesh = feflow_model.mesh # Print information about the mesh. print(mesh) .. rst-class:: sphx-glr-script-out .. code-block:: none UnstructuredGrid (0x7b26fcde4160) N Cells: 6260 N Points: 3228 X Bounds: -2.500e+02, 9.000e+02 Y Bounds: 3.500e+02, 6.500e+02 Z Bounds: 0.000e+00, 0.000e+00 N Arrays: 34 .. GENERATED FROM PYTHON SOURCE LINES 32-33 2. Plot the temperature simulated in FEFLOW on the mesh. .. GENERATED FROM PYTHON SOURCE LINES 33-34 .. code-block:: Python fig = ot.plot.contourf(mesh, "P_TEMP", show_edges=True) .. image-sg:: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_001.png :alt: plot B feflowlib BC mesh :srcset: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_001.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 35-39 3. As the FEFLOW data now are a pyvista.UnstructuredGrid, all pyvista functionalities can be applied to it. Further information can be found at https://docs.pyvista.org/version/stable/user-guide/simple.html. For example it can be saved as a VTK Unstructured Grid File (\*.vtu). This allows to use the FEFLOW model for ``OGS`` simulation or to observe it in ``Paraview```. .. GENERATED FROM PYTHON SOURCE LINES 39-40 .. code-block:: Python pv.save_meshio(temp_dir / "HT_mesh.vtu", mesh) .. GENERATED FROM PYTHON SOURCE LINES 41-42 4. Run the FEFLOW model in OGS. .. GENERATED FROM PYTHON SOURCE LINES 42-45 .. code-block:: Python feflow_model.setup_prj(end_time=1e11, time_stepping=[(1, 1e10)]) feflow_model.run() ms = ot.MeshSeries(temp_dir / "HT_model.pvd") .. rst-class:: sphx-glr-script-out .. code-block:: none OGS finished with project file /tmp/tmpeim6fyw1converted_models/HT_model.prj. Execution took 1.472346544265747 s Project file written to output. .. GENERATED FROM PYTHON SOURCE LINES 46-47 5. Plot the temperature simulated in OGS. .. GENERATED FROM PYTHON SOURCE LINES 47-49 .. code-block:: Python ogs_sim_res = ms.mesh(ms.timesteps[-1]) ot.plot.contourf(ogs_sim_res, ot.variables.temperature) .. image-sg:: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_002.png :alt: plot B feflowlib BC mesh :srcset: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_002.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none
.. GENERATED FROM PYTHON SOURCE LINES 50-53 6. The boundary meshes are manipulated to alter the model. The original boundary conditions are shown in this example: :ref:`sphx_glr_auto_examples_howto_conversions_plot_F_feflowlib_HT_simulation.py` 6.1 The Dirichlet boundary conditions for the hydraulic head are set to 0. Therefore, no water flows from the left to the right edge. .. GENERATED FROM PYTHON SOURCE LINES 53-55 .. code-block:: Python bc_flow = feflow_model.subdomains["P_BC_FLOW"]["P_BC_FLOW"] feflow_model.subdomains["P_BC_FLOW"]["P_BC_FLOW"][bc_flow == 10] = 0 .. GENERATED FROM PYTHON SOURCE LINES 56-57 6.2 Overwrite the new boundary conditions and run the model. .. GENERATED FROM PYTHON SOURCE LINES 57-58 .. code-block:: Python feflow_model.run(overwrite=True) .. rst-class:: sphx-glr-script-out .. code-block:: none OGS finished with project file /tmp/tmpeim6fyw1converted_models/HT_model.prj. Execution took 1.4022948741912842 s Project file written to output. .. GENERATED FROM PYTHON SOURCE LINES 59-60 6.3 The corresponding simulation results look like. .. GENERATED FROM PYTHON SOURCE LINES 60-64 .. code-block:: Python ms = ot.MeshSeries(temp_dir / "HT_model.pvd") # Read the last timestep: ogs_sim_res = ms.mesh(ms.timesteps[-1]) ot.plot.contourf(ogs_sim_res, ot.variables.temperature) .. image-sg:: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_003.png :alt: plot B feflowlib BC mesh :srcset: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_003.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none
.. GENERATED FROM PYTHON SOURCE LINES 65-66 6.4 Create a new boundary mesh and overwrite the existing subdomains with this boundary mesh. .. GENERATED FROM PYTHON SOURCE LINES 66-78 .. code-block:: Python assign_bulk_ids(mesh) # Get the points of the bulk mesh to build a new boundary mesh. wanted_pts = [1492, 1482, 1481, 1491, 1479, 1480, 1839, 1840] new_bc = mesh.extract_points( [pt in wanted_pts for pt in mesh["bulk_node_ids"]], adjacent_cells=False, include_cells=False, ) new_bc["bulk_node_ids"] = np.array(wanted_pts, dtype=np.uint64) # Define the temperature values of these points. new_bc["P_BC_HEAT"] = np.array([300] * len(wanted_pts), dtype=np.float64) feflow_model.subdomains["P_BC_HEAT"] = new_bc .. GENERATED FROM PYTHON SOURCE LINES 79-80 7. Run the new model and plot the simulation results. .. GENERATED FROM PYTHON SOURCE LINES 80-84 .. code-block:: Python feflow_model.run(overwrite=True) ms = ot.MeshSeries(temp_dir / "HT_model.pvd") ogs_sim_res = ms.mesh(ms.timesteps[-1]) ot.plot.contourf(ogs_sim_res, ot.variables.temperature) .. image-sg:: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_004.png :alt: plot B feflowlib BC mesh :srcset: /auto_examples/howto_conversions/images/sphx_glr_plot_B_feflowlib_BC_mesh_004.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none OGS finished with project file /tmp/tmpeim6fyw1converted_models/HT_model.prj. Execution took 1.419546365737915 s Project file written to output.
.. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 5.955 seconds) .. _sphx_glr_download_auto_examples_howto_conversions_plot_B_feflowlib_BC_mesh.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_B_feflowlib_BC_mesh.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_B_feflowlib_BC_mesh.py ` .. container:: sphx-glr-download sphx-glr-download-zip :download:`Download zipped: plot_B_feflowlib_BC_mesh.zip `