Contents
ogs5py Quickstart

ogs5py is A python-API for the OpenGeoSys 5 scientific modeling package.
Installation
The package can be installed via pip. On Windows you can install WinPython to get Python and pip running.
pip install ogs5py
Or with conda:
conda install ogs5py
Citation
If you are using ogs5py in your publication please cite our paper:
Müller, S., Zech, A. and Heße, F.: ogs5py: A Python-API for the OpenGeoSys 5 Scientific Modeling Package. Groundwater, 59: 117-122. https://doi.org/10.1111/gwat.13017, 2021.
You can cite the Zenodo code publication of ogs5py by:
Sebastian Müller. GeoStat-Framework/ogs5py. Zenodo. https://doi.org/10.5281/zenodo.2546767
If you want to cite a specific version, have a look at the Zenodo site.
Further Information
General homepage: https://www.opengeosys.org/ogs-5
OGS5 Repository: https://github.com/ufz/ogs5
Keyword documentation: https://ogs5-keywords.netlify.com
OGS5 Benchmarks: https://github.com/ufz/ogs5-benchmarks
ogs5py Benchmarks: https://github.com/GeoStat-Framework/ogs5py_benchmarks
Pumping Test Example
In the following a simple transient pumping test is simulated on a radial symmetric mesh. The point output at the observation well is plotted afterwards.
from ogs5py import OGS, specialrange, generate_time
from matplotlib import pyplot as plt
# discretization and parameters
time = specialrange(0, 3600, 50, typ="cub")
rad = specialrange(0, 1000, 100, typ="cub")
obs = rad[21]
angles = 32
storage = 1e-3
transmissivity = 1e-4
rate = -1e-3
# model setup
model = OGS(task_root="pump_test", task_id="model")
# generate a radial mesh and geometry ("boundary" polyline)
model.msh.generate("radial", dim=2, rad=rad, angles=angles)
model.gli.generate("radial", dim=2, rad_out=rad[-1], angles=angles)
model.gli.add_points([0.0, 0.0, 0.0], "pwell")
model.gli.add_points([obs, 0.0, 0.0], "owell")
model.bc.add_block( # boundary condition
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["POLYLINE", "boundary"],
DIS_TYPE=["CONSTANT", 0.0],
)
model.st.add_block( # source term
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["POINT", "pwell"],
DIS_TYPE=["CONSTANT_NEUMANN", rate],
)
model.ic.add_block( # initial condition
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE="DOMAIN",
DIS_TYPE=["CONSTANT", 0.0],
)
model.mmp.add_block( # medium properties
GEOMETRY_DIMENSION=2,
STORAGE=[1, storage],
PERMEABILITY_TENSOR=["ISOTROPIC", transmissivity],
)
model.num.add_block( # numerical solver
PCS_TYPE="GROUNDWATER_FLOW",
LINEAR_SOLVER=[2, 5, 1e-14, 1000, 1.0, 100, 4],
)
model.out.add_block( # point observation
PCS_TYPE="GROUNDWATER_FLOW",
NOD_VALUES="HEAD",
GEO_TYPE=["POINT", "owell"],
DAT_TYPE="TECPLOT",
)
model.pcs.add_block( # set the process type
PCS_TYPE="GROUNDWATER_FLOW", NUM_TYPE="NEW"
)
model.tim.add_block( # set the timesteps
PCS_TYPE="GROUNDWATER_FLOW",
**generate_time(time)
)
model.write_input()
model.run_model()

OGS5 executable
To obtain an OGS5 executable, ogs5py
brings a download routine download_ogs
:
from ogs5py import download_ogs
download_ogs()
Then a executable is stored in the ogs5py config path and will be called when a model is run.
You can pass a version
statement to the download_ogs
routine, to
obtain a specific version (5.7, 5.7.1 (win only) and 5.8).
For OGS 5.7 there are executables for Windows/Linux and MacOS.
For “5.8” there are no MacOS pre-builds.
If you have compiled your own OGS5 version, you can add your executable
to the ogs5py config path with add_exe
:
from ogs5py import add_exe
add_exe("path/to/your/ogs/exe")
Otherwise you need to specify the path to the executable within the run command:
model.run_model(ogs_exe="path/to/ogs")
Requirements
License
ogs5py Tutorials
In the following you will find several Tutorials on how to use ogs5py to explore its whole beauty and power.
Tutorial 1: A pumping test
This is a minimal example on how to setup a pumping test with ogs5py. The result was plotted against the analytical solution.
In this example we use the generate_time
function, to use
an array of time points for the time stepping definition.
model.tim.add_block(**generate_time(time))
is equivalent to:
model.tim.add_block(
TIME_START=0,
TIME_END=time[-1],
TIME_STEPS=[
[1, time[0]],
[1, time[1]],
[1, time[2]],
# ...
],
)
The script:
import anaflow as ana
from ogs5py import OGS, specialrange, generate_time
from matplotlib import pyplot as plt
# discretization and parameters
time = specialrange(0, 3600, 50, typ="cub")
rad = specialrange(0, 1000, 100, typ="cub")
obs = rad[21]
angles = 32
storage = 1e-3
transmissivity = 1e-4
rate = -1e-3
# model setup
model = OGS(task_root="pump_test", task_id="model")
model.pcs.add_block( # set the process type
PCS_TYPE="GROUNDWATER_FLOW", NUM_TYPE="NEW"
)
# generate a radial mesh and geometry ("boundary" polyline)
model.msh.generate("radial", dim=2, rad=rad, angles=angles)
model.gli.generate("radial", dim=2, rad_out=rad[-1], angles=angles)
model.gli.add_points([0.0, 0.0, 0.0], "pwell")
model.gli.add_points([obs, 0.0, 0.0], "owell")
model.bc.add_block( # boundary condition
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["POLYLINE", "boundary"],
DIS_TYPE=["CONSTANT", 0.0],
)
model.ic.add_block( # initial condition
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE="DOMAIN",
DIS_TYPE=["CONSTANT", 0.0],
)
model.st.add_block( # source term
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["POINT", "pwell"],
DIS_TYPE=["CONSTANT_NEUMANN", rate],
)
model.mmp.add_block( # medium properties
GEOMETRY_DIMENSION=2,
STORAGE=[1, storage],
PERMEABILITY_TENSOR=["ISOTROPIC", transmissivity],
)
model.num.add_block( # numerical solver
PCS_TYPE="GROUNDWATER_FLOW",
LINEAR_SOLVER=[2, 5, 1e-14, 1000, 1.0, 100, 4]
)
model.out.add_block( # point observation
PCS_TYPE="GROUNDWATER_FLOW",
NOD_VALUES="HEAD",
GEO_TYPE=["POINT", "owell"],
DAT_TYPE="TECPLOT",
)
model.tim.add_block( # set the timesteps
PCS_TYPE="GROUNDWATER_FLOW",
**generate_time(time) # generate input from time-series
)
model.write_input()
success = model.run_model()
print("success:", success)
# observation
point = model.readtec_point(pcs="GROUNDWATER_FLOW")
time = point["owell"]["TIME"]
head = point["owell"]["HEAD"]
# analytical solution
head_ana = ana.theis(time, obs, storage, transmissivity, rate=rate)
# comparisson plot
plt.scatter(time, head, color="k", label="simulated, r={:04.2f}m".format(obs))
plt.plot(time, head_ana, label="analytical solution")
plt.xscale("symlog", linthreshx=10, subsx=range(1, 10))
plt.xlim([0, 1.1 * time[-1]])
plt.xlabel("time in s")
plt.ylabel("head in m")
plt.legend()
plt.show()
# show mesh
model.msh.show()


Tutorial 2: Interaction with GSTools
In this example we are generating a log-normal dirtributed conductivity field on a generated mesh with the aid of GSTools. and perform a steady pumping test.
import numpy as np
from ogs5py import OGS, by_id, show_vtk
from gstools import SRF, Gaussian
# covariance model for conductivity field
cov_model = Gaussian(dim=3, var=2, len_scale=10, anis=[1, 0.2])
srf = SRF(model=cov_model, mean=-9, seed=1000)
# model setup
model = OGS(task_root="test_het_3D", task_id="model", output_dir="out")
model.pcs.add_block( # set the process type
PCS_TYPE="GROUNDWATER_FLOW", NUM_TYPE="NEW", TIM_TYPE="STEADY"
)
# generate a radial 3D mesh and conductivity field
model.msh.generate(
"radial", dim=3, angles=64, rad=np.arange(101), z_arr=-np.arange(11)
)
cond = np.exp(srf.mesh(model.msh))
model.mpd.add(name="conductivity")
model.mpd.add_block( # edit recent mpd file
MSH_TYPE="GROUNDWATER_FLOW",
MMP_TYPE="PERMEABILITY",
DIS_TYPE="ELEMENT",
DATA=by_id(cond),
)
model.mmp.add_block( # permeability, storage and porosity
GEOMETRY_DIMENSION=3, PERMEABILITY_DISTRIBUTION=model.mpd.file_name
)
model.gli.generate("radial", dim=3, angles=64, rad_out=100, z_size=-10)
model.gli.add_polyline("pwell", [[0, 0, 0], [0, 0, -10]])
for srf in model.gli.SURFACE_NAMES: # set boundary condition
model.bc.add_block(
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["SURFACE", srf],
DIS_TYPE=["CONSTANT", 0.0],
)
model.st.add_block( # set pumping condition at the pumpingwell
PCS_TYPE="GROUNDWATER_FLOW",
PRIMARY_VARIABLE="HEAD",
GEO_TYPE=["POLYLINE", "pwell"],
DIS_TYPE=["CONSTANT_NEUMANN", 1.0e-3],
)
model.num.add_block( # numerical solver
PCS_TYPE="GROUNDWATER_FLOW",
LINEAR_SOLVER=[2, 5, 1.0e-14, 1000, 1.0, 100, 4],
)
model.out.add_block( # set the outputformat
PCS_TYPE="GROUNDWATER_FLOW",
NOD_VALUES="HEAD",
GEO_TYPE="DOMAIN",
DAT_TYPE="VTK",
)
model.write_input()
success = model.run_model()
model.msh.show(show_cell_data={"Conductivity": cond}, log_scale=True)
files = model.output_files(pcs="GROUNDWATER_FLOW", typ="VTK")
show_vtk(files[-1], log_scale=True) # show the last time-step


Tutorial 2: Interaction with pygmsh
In this example we are generating different meshes with the aid of pygmsh and gmsh
import numpy as np
import pygmsh
from ogs5py import OGS
with pygmsh.geo.Geometry() as geom:
poly = geom.add_polygon(
[
[+0.0, +0.5],
[-0.1, +0.1],
[-0.5, +0.0],
[-0.1, -0.1],
[+0.0, -0.5],
[+0.1, -0.1],
[+0.5, +0.0],
[+0.1, +0.1],
],
mesh_size=0.05,
)
geom.twist(
poly,
translation_axis=[0, 0, 1],
rotation_axis=[0, 0, 1],
point_on_axis=[0, 0, 0],
angle=np.pi / 3,
)
mesh = geom.generate_mesh()
model = OGS()
# generate example above
model.msh.import_mesh(mesh, import_dim=3)
model.msh.show()
# generate a predefined grid adapter in 2D
model.msh.generate("grid_adapter2D", in_mat=1, out_mat=0, fill=True)
model.msh.show(show_material_id=True)
# generate a predefined grid adapter in 3D
model.msh.generate("grid_adapter3D", in_mat=1, out_mat=0, fill=True)
model.msh.show(show_material_id=True)
# generate a predefined block adapter in 3D
model.msh.generate("block_adapter3D", xy_dim=5.0, z_dim=1.0, in_res=1)
model.msh.show(show_element_id=True)




ogs5py API
Purpose
ogs5py is A python-API for the OpenGeoSys 5 scientific modeling package.
The following functionalities are directly provided on module-level.
Subpackages
ogs5py subpackage providing the file classes. |
|
ogs5py subpackage providing reader for the ogs5 output. |
|
ogs5py subpackage providing tools. |
Classes
OGS model Base Class
Class to setup an ogs model
|
Class for an OGS5 model. |
File Classes
Classes for all OGS5 Files. See: ogs5py.fileclasses
|
Class for the ogs ASC file. |
|
Class for the ogs BOUNDARY CONDITION file. |
|
Class for the ogs COMMUNICATION TABLE file. |
|
Class for the ogs MPI DOMAIN DECOMPOSITION file. |
|
Class for the ogs FUNCTION file. |
|
Class for the ogs GEOCHEMICAL THERMODYNAMIC MODELING COUPLING file. |
|
Class for GEMS3K input file. |
|
Class for the ogs GEOMETRY file. |
|
Class for an external definition for the ogs GEOMETRY file. |
|
Class for the ogs INITIAL_CONDITION file. |
|
Class for the ogs RESTART file, if the DIS_TYPE in IC is set to RESTART. |
|
Class for the ogs KINETRIC REACTION file. |
|
Class for the ogs COMPONENT_PROPERTIES file. |
|
Class for the ogs FLUID PROPERTY file. |
|
Class for the ogs MEDIUM_PROPERTIES file. |
|
Class for the ogs MEDIUM_PROPERTIES_DISTRIBUTED file. |
|
Class for a multi layer mesh file that contains multiple '#FEM_MSH' Blocks. |
|
Class for the ogs SOLID_PROPERTIES file. |
|
Class for the ogs NUMERICS file. |
|
Class for the ogs OUTPUT file. |
|
Class for the ogs PROCESS file. |
|
Class for the ogs Particle file, if the PCS TYPE is RANDOM_WALK. |
|
Class for the ogs PHREEQC interface file. |
|
Class for the ogs PHREEQC dat file. |
|
Class for the ogs REACTION_INTERFACE file. |
|
Class for the ogs USER DEFINED TIME CURVES file. |
|
Class for the ogs SOURCE_TERM file. |
|
Class for the ogs TIME_STEPPING file. |
Functions
Geometric
Geometric routines
|
Providing a transformation function to deform a given mesh. |
Searching
Routine to search for a valid ogs id in a directory
|
Search for OGS model names in the given path. |
Formatting
Routines to format/generate data in the right way for the input
|
Return a flattend array side-by-side with the array-element ids. |
|
Calculation of special point ranges. |
|
Return a dictionary for the ".tim" file. |
Downloading
Routine to download OGS5.
|
Download the OGS5 executable. |
|
Add an OGS5 exe to |
Reset all downloads in |
|
Standard config path for ogs5py. |
Plotting
Routine to download OGS5.
|
Display a given mesh colored by its material ID. |
Information
all ogs file extensions |
|
PCS types |
|
primary variables per PCS |
Changelog
All notable changes to ogs5py will be documented in this file.
1.3.0 - 2023-04
See #18
Enhancements
move to
src/
base package structureuse hatchling as build backend
drop py36 support
added archive support
simplify documentation
Bugfixes
remove usage of deprecated
np.int
1.2.2 - 2022-05-25
Bugfixes
1.2.1 - 2022-05-15
Enhancements
MSH.import_mesh
can handlemeshio.Mesh
as input now #13
Changes
pygmsh
support was removed. You can’t use pygmsh Geometry objects to generate meshes anymore. Please generate beforehand and import the generated mesh. Other generators are usinggmsh
directly now. #13
1.2.0 - 2022-05-15
Enhancements
move to a
pyproject.toml
based installation: d5ea756move from
develop
/master
branches to a singlemain
branchuse GitHub Actions for CI: b6811ce
use f-strings where possible #11
simplified documentation #11
added changelog to documentation #11
added citation file and paper reference #11
use Python 3 style classes #11
Changes
downlaod_ogs
only supports version “5.7”, “5.7.1” and “5.8” since the CI for OGS5 was shut down: 8b1cc91
Bugfixes
1.1.1 - 2020-04-02
Bugfixes
check if
__version__
is present (only if installed)
1.1.0 - 2020-03-22
Bugfixes
meshio 4 was not compatible
fixed integer type in exporting meshes with element/material IDs
better check for OGS5 success on Windows
Changes
drop py2.7 support
1.0.5 - 2019-11-18
Bugfixes
Additions
1.0.4 - 2019-09-10
Bugfixes
Additions
new sub-keywords for OUT (added to OGS5 in Aug 19) when using TECPLOT (TECPLOT_ELEMENT_OUTPUT_CELL_CENTERED, TECPLOT_ZONES_FOR_MG): ebcb22a
Changes
RFR Class was refactored to allow multiple variables: 3c1c445
1.0.3 - 2019-08-23
Bugfixes
MSH.show
TempFile was not working on Windows: c0d0960
1.0.2 - 2019-08-22
Bugfixes
1.0.1 - 2019-08-22
Bugfixes
download_ogs(version="latest" build="PETSC")
was not working: 552503b
1.0.0 - 2019-08-22
Bugfixes
GLI.add_polyline
now allows integer coordinates for points: bf5d684MSH.centroids
are now calculated as center of mass instead of center of element nodes: b0708a6MSH.show
was not working: 6a0489bOGS.run_model
has now a better check for OGS success: 143d0abGMSH
interface was updated to new meshio-API: d3e0594RFR
file was not written: 41e55f3BC
new sub-key TIME_INTERVAL was missing: 94ec5c5
Additions
download_ogs
downloads a system dependent OGS5 executable: ede32e4add_exe
add a self compiled OGS5 executable: ede32e4MSH.import_mesh
now allows the import of material_id and element_id if given as cell_data in the external mesh: 00a77faMSH.export_mesh
now automatically exports material_id (already the case before) and element_id. Also you can now export additionalpoint_data
andfield_data
: 00a77faNew method
MSH.set_material_id
to set the material IDs for specific elements: 4b11c6aMSH.show
now can show additional cell_data: ffd7604New routine
show_vtk
to show vtk output with mayavi: f640c19New method
OGS.output_files
to get a list of output files: 2f5f102New attribute
file_name
for files: 632c2e7BlockFile: new method
append_to_block
: efc9aacOGS.gen_script
now allowes multiple subkeys: 2cd344b
Changes
MSH.export_mesh
argumentadd_data_by_id
renamed tocell_data_by_id
: 00a77faOGS.run_model
argumentogs_root
renamed toogs_exe
: 6fcdb61Files that can occure multiple times (mpd, rfr, …) are better handled now: 4a9c9d2
ogs5py
is now licensed under the MIT license: ae96c0eExtra named files now get their name by keyword
name
: 632c2e7
0.6.5 - 2019-07-05
Bugfixes
gli.add_polyline: Adding polyline by given point IDs was not possible: 3ec23af
Additions
New swap_axis routine in msh and gli: You can now easily swap axis of a mesh. If you have generated a 2D mesh in x-y you can get a x-z cross-section by swapping the y and z axis: 3ec23af
0.6.4 - 2019-05-16
Bugfixes
Additions
0.6.3 - 2019-03-21
Bugfixes
The used method os.makedirs has no keyword argument ‘exist_ok’ in python 2.7, so we prevent using it. See: 40fea36
0.6.2 - 2019-03-21
Bugfixes
The vtk reading routine could not read multiple scalar cell data. See: 568e7be
0.6.1 - 2019-01-22
Bugfixes
The BlockFile reading routine was cutting of the given keys. See: d82fd30
0.6.0 - 2019-01-22
First release of ogs5py.