.. DO NOT EDIT. .. THIS FILE WAS AUTOMATICALLY GENERATED BY SPHINX-GALLERY. .. TO MAKE CHANGES, EDIT THE SOURCE PYTHON FILE: .. "auto_examples/howto_feflowlib/plot_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_simulation.py: How to convert a FEFLOW model and simulate it in OGS. ===================================================== .. sectionauthor:: Julian Heinze (Helmholtz Centre for Environmental Research GmbH - UFZ) In this example we show how a simple 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-33 .. code-block:: default import tempfile import xml.etree.ElementTree as ET from pathlib import Path import ifm_contrib as ifm from ogs6py import ogs from pyvista import read from ogstools.feflowlib import ( convert_properties_mesh, extract_cell_boundary_conditions, setup_prj_file, steady_state_diffusion, ) from ogstools.feflowlib.examples import path_box_Neumann from ogstools.feflowlib.tools import ( extract_point_boundary_conditions, get_material_properties, ) .. GENERATED FROM PYTHON SOURCE LINES 34-35 1. Load a FEFLOW model (.fem) as a FEFLOW document and convert it. .. GENERATED FROM PYTHON SOURCE LINES 35-44 .. code-block:: default feflow_model = ifm.loadDocument(path_box_Neumann) pyvista_mesh = convert_properties_mesh(feflow_model) pyvista_mesh.plot( show_edges=True, off_screen=True, scalars="P_HEAD", cpos=[0, 1, 0.5] ) print(pyvista_mesh) path_writing = Path(tempfile.mkdtemp("feflow_test_simulation")) path_mesh = path_writing / "boxNeumann.vtu" pyvista_mesh.save(str(path_mesh)) .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_001.png :alt: plot simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_001.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-script-out .. code-block:: none UnstructuredGrid (0x7f59355e1880) 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 45-46 2. Extract the point and cell boundary conditions and write them to a temporary directory. .. GENERATED FROM PYTHON SOURCE LINES 46-58 .. code-block:: default point_BC_dict = extract_point_boundary_conditions(path_writing, pyvista_mesh) # Since there can be multiple point boundary conditions on the bulk mesh, # they are saved and plotted iteratively. for path, boundary_condition in point_BC_dict.items(): boundary_condition.save(path) boundary_condition.plot() path_topsurface, topsurface = extract_cell_boundary_conditions( path_mesh, 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) .. rst-class:: sphx-glr-horizontal * .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_002.png :alt: plot simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_002.png :class: sphx-glr-multi-img * .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_003.png :alt: plot simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_003.png :class: sphx-glr-multi-img .. GENERATED FROM PYTHON SOURCE LINES 59-60 3. Setup a prj-file to run a OGS-simulation .. GENERATED FROM PYTHON SOURCE LINES 60-81 .. code-block:: default path_prjfile = str(path_mesh.with_suffix(".prj")) prjfile = ogs.OGS(PROJECT_FILE=str(path_prjfile)) # Get the template prj-file configurations for a steady state diffusion process model = steady_state_diffusion( path_writing / "sim_boxNeumann", prjfile, ) # Include the mesh specific configurations to the template. model = setup_prj_file( path_mesh, pyvista_mesh, get_material_properties(pyvista_mesh, "P_CONDX"), "steady state diffusion", model, ) # The model must be written before it can be run. model.write_input(str(path_prjfile)) # Simply print the prj-file as an example. read_model = ET.parse(str(path_prjfile)) root = read_model.getroot() ET.dump(root) .. 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/tmpiepoz_brfeflow_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 P_IOFLOW MeshElement topsurface_boxNeumann P_IOFLOW P_SOUF MeshElement topsurface_boxNeumann P_SOUF HEAD_OGS 1 1 p0 Dirichlet P_BC_FLOW P_BC_FLOW Neumann P_BCFLOW_2ND P_BCFLOW_2ND Neumann topsurface_boxNeumann P_IOFLOW Volumetric topsurface_boxNeumann P_SOUF 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 82-83 4. Run the model .. GENERATED FROM PYTHON SOURCE LINES 83-84 .. code-block:: default model.run_model(logfile=str(path_writing / "out.log")) .. rst-class:: sphx-glr-script-out .. code-block:: none OGS finished with project file /tmp/tmpiepoz_brfeflow_test_simulation/boxNeumann.prj. Execution took 0.29024648666381836 s .. GENERATED FROM PYTHON SOURCE LINES 85-86 5. Read the results and plot them. .. GENERATED FROM PYTHON SOURCE LINES 86-91 .. code-block:: default ogs_sim_res = read(str(path_writing / "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] ) .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_004.png :alt: plot simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_004.png :class: sphx-glr-single-img .. GENERATED FROM PYTHON SOURCE LINES 92-93 6. Calculate the difference to the FEFLOW simulation and plot it. .. GENERATED FROM PYTHON SOURCE LINES 93-98 .. code-block:: default 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] ) .. image-sg:: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_005.png :alt: plot simulation :srcset: /auto_examples/howto_feflowlib/images/sphx_glr_plot_simulation_005.png :class: sphx-glr-single-img .. rst-class:: sphx-glr-timing **Total running time of the script:** (0 minutes 1.269 seconds) .. _sphx_glr_download_auto_examples_howto_feflowlib_plot_simulation.py: .. only:: html .. container:: sphx-glr-footer sphx-glr-footer-example .. container:: sphx-glr-download sphx-glr-download-python :download:`Download Python source code: plot_simulation.py ` .. container:: sphx-glr-download sphx-glr-download-jupyter :download:`Download Jupyter notebook: plot_simulation.ipynb `