.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/howto_feflowlib/plot_H_simulation.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_feflowlib_plot_H_simulation.py: Hydraulic model - conversion and simulation =========================================== .. sectionauthor:: Julian Heinze (Helmholtz Centre for Environmental Research GmbH - UFZ) In this example we show how a simple flow/hydraulic FEFLOW model can be converted to a pyvista.UnstructuredGrid and then be simulated in OGS. .. GENERATED FROM PYTHON SOURCE LINES 12-13 0. Necessary imports .. GENERATED FROM PYTHON SOURCE LINES 13-34 .. code-block:: Python import tempfile import xml.etree.ElementTree as ET from pathlib import Path import ifm_contrib as ifm import pyvista as pv from ogs6py import ogs from ogstools.examples import feflow_model_box_Neumann from ogstools.feflowlib import ( convert_properties_mesh, extract_cell_boundary_conditions, setup_prj_file, steady_state_diffusion, ) from ogstools.feflowlib.tools import ( extract_point_boundary_conditions, get_material_properties, ) from ogstools.meshlib import MeshSeries .. GENERATED FROM PYTHON SOURCE LINES 35-37 1. Load a FEFLOW model (.fem) as a FEFLOW document, convert and save it. More details on how the conversion function works can be found here: :py:mod:`ogstools.feflowlib.convert_properties_mesh`. .. GENERATED FROM PYTHON SOURCE LINES 37-52 .. code-block:: Python feflow_model = ifm.loadDocument(str(feflow_model_box_Neumann)) pyvista_mesh = convert_properties_mesh(feflow_model) pv.global_theme.colorbar_orientation = "vertical" pyvista_mesh.plot( show_edges=True, off_screen=True, scalars="P_HEAD", cpos=[0, 1, 0.5], scalar_bar_args={"position_x": 0.1, "position_y": 0.25}, ) print(pyvista_mesh) temp_dir = Path(tempfile.mkdtemp("feflow_test_simulation")) feflow_mesh_file = temp_dir / "boxNeumann.vtu" pyvista_mesh.save(feflow_mesh_file) .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_001.png :alt: plot H simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_001.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /builds/ogs/tools/ogstools/docs/auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_001.vtksz .. rst-class:: sphx-glr-script-out .. code-block:: none UnstructuredGrid (0x7fb14ee61e80) N Cells: 11462 N Points: 6768 X Bounds: 0.000e+00, 1.000e+02 Y Bounds: 0.000e+00, 1.000e+02 Z Bounds: -1.000e+02, 0.000e+00 N Arrays: 23 .. GENERATED FROM PYTHON SOURCE LINES 53-54 2. Extract the point conditions (see: :py:mod:`ogstools.feflowlib.extract_point_boundary_conditions`). .. GENERATED FROM PYTHON SOURCE LINES 54-69 .. code-block:: Python point_BC_dict = extract_point_boundary_conditions(temp_dir, pyvista_mesh) # Since there can be multiple point boundary conditions on the bulk mesh, # they are saved and plotted iteratively. plotter = pv.Plotter(shape=(len(point_BC_dict), 1)) for i, (path, boundary_condition) in enumerate(point_BC_dict.items()): boundary_condition.save(path) plotter.subplot(i, 0) plotter.add_mesh(boundary_condition, scalars=Path(path).stem) plotter.show() path_topsurface, topsurface = extract_cell_boundary_conditions( feflow_mesh_file, pyvista_mesh ) # On the topsurface can be cell based boundary condition. # The boundary conditions on the topsurface of the model are required for generalization. topsurface.save(path_topsurface) .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_002.png :alt: plot H simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_002.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /builds/ogs/tools/ogstools/docs/auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_002.vtksz .. GENERATED FROM PYTHON SOURCE LINES 70-71 3. Setup a prj-file (see: :py:mod:`ogstools.feflowlib.setup_prj_file`) to run a OGS-simulation. .. GENERATED FROM PYTHON SOURCE LINES 71-91 .. code-block:: Python path_prjfile = feflow_mesh_file.with_suffix(".prj") prjfile = ogs.OGS(PROJECT_FILE=path_prjfile) # Get the template prj-file configurations for a steady state diffusion process ssd_model = steady_state_diffusion( temp_dir / "sim_boxNeumann", prjfile, ) # Include the mesh specific configurations to the template. model = setup_prj_file( bulk_mesh_path=feflow_mesh_file, mesh=pyvista_mesh, material_properties=get_material_properties(pyvista_mesh, "P_CONDX"), process="steady state diffusion", model=ssd_model, ) # The model must be written before it can be run. model.write_input(path_prjfile) # Simply print the prj-file as an example. model_prjfile = ET.parse(path_prjfile) ET.dump(model_prjfile) .. rst-class:: sphx-glr-script-out .. code-block:: none boxNeumann.vtu topsurface_boxNeumann.vtu P_BC_FLOW.vtu P_BCFLOW_2ND.vtu SteadyStateDiffusion STEADY_STATE_DIFFUSION 2 HEAD_OGS diffusion Constant 1.1574074074074073e-05 reference_temperature Constant 293.15 basic_picard DeltaX NORM2 1e-15 BackwardEuler SingleStep VTK /tmp/tmpx43eawxyfeflow_test_simulation/sim_boxNeumann 1 1 p0 Constant 0 P_BC_FLOW MeshNode P_BC_FLOW P_BC_FLOW P_BCFLOW_2ND MeshNode P_BCFLOW_2ND P_BCFLOW_2ND HEAD_OGS 1 1 p0 Dirichlet P_BC_FLOW P_BC_FLOW Neumann P_BCFLOW_2ND P_BCFLOW_2ND basic_picard Picard 10 general_linear_solver general_linear_solver -i cg -p jacobi -tol 1e-6 -maxiter 100000 CG DIAGONAL 100000 1e-6 sd -sd_ksp_type cg -sd_pc_type bjacobi -sd_ksp_rtol 1e-16 -sd_ksp_max_it 10000 .. GENERATED FROM PYTHON SOURCE LINES 92-93 4. Run the model .. GENERATED FROM PYTHON SOURCE LINES 93-94 .. code-block:: Python model.run_model(logfile=temp_dir / "out.log") .. rst-class:: sphx-glr-script-out .. code-block:: none OGS finished with project file /tmp/tmpx43eawxyfeflow_test_simulation/boxNeumann.prj. Execution took 0.41233348846435547 s .. GENERATED FROM PYTHON SOURCE LINES 95-96 5. Read the results and plot them. .. GENERATED FROM PYTHON SOURCE LINES 96-111 .. code-block:: Python ms = MeshSeries(temp_dir / "sim_boxNeumann.pvd") # Read the last timestep: ogs_sim_res = ms.read(ms.timesteps[-1]) """ It is also possible to read the file directly with pyvista: ogs_sim_res = pv.read(temp_dir / "sim_boxNeumann_ts_1_t_1.000000.vtu") """ ogs_sim_res.plot( show_edges=True, off_screen=True, scalars="HEAD_OGS", cpos=[0, 1, 0.5], scalar_bar_args={"position_x": 0.1, "position_y": 0.25}, ) .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_003.png :alt: plot H simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_003.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /builds/ogs/tools/ogstools/docs/auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_003.vtksz .. GENERATED FROM PYTHON SOURCE LINES 112-113 6. Calculate the difference to the FEFLOW simulation and plot it. .. GENERATED FROM PYTHON SOURCE LINES 113-122 .. code-block:: Python diff = pyvista_mesh["P_HEAD"] - ogs_sim_res["HEAD_OGS"] pyvista_mesh["diff_HEAD"] = diff pyvista_mesh.plot( show_edges=True, off_screen=True, scalars="diff_HEAD", cpos=[0, 1, 0.5], scalar_bar_args={"position_x": 0.1, "position_y": 0.25}, ) .. tab-set:: .. tab-item:: Static Scene .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_004.png :alt: plot H simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_004.png :class: sphx-glr-single-img .. tab-item:: Interactive Scene .. offlineviewer:: /builds/ogs/tools/ogstools/docs/auto_examples/howto_feflowlib/images/sphx_glr_plot_H_simulation_004.vtksz .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 2.160 seconds) .. _sphx_glr_download_auto_examples_howto_feflowlib_plot_H_simulation.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_H_simulation.ipynb ` .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_H_simulation.py `