HyVR package

hyvr.grid module

Grid class

A module containing some classes and function useful to work with structured grids.
Usage:

Import as a normal python module.

Version:

0.1 , 01-09-2016 : Forked from Alessandro Comunian

Authors:

Jeremy P. Bennett

Notes:

The grids for the moment are always considered as 3D.

class hyvr.grid.Grid(ox=0.0, oy=0.0, oz=0.0, dx=1.0, dy=1.0, dz=1.0, nx=200, ny=200, nz=10, gtype='points', gname='image', periodicity=False)[source]

Bases: object

Grid class

A simple class that contains the Origin, Delta and Size of a simulation. It can also be used as container for some information contained in a VTK structured grid file.

Notes

  • By default, the size of the grid is considered as points.

See also

vtknumpy

cart_coords()[source]

Get x,y,z coordinates in a tuple

Parameters:self – An instance of the Grid class
Returns:A Tuple containing the x,y and z-coordinates of the grid
cart_coords2d()[source]

Get x,y coordinates in a tuple

Parameters:self – An instance of the Grid class
Returns:A tuple containing the x,y-coordinates of the grid
cells

Number of cells

cellsize_2d()[source]

Compute the cell size in 2D

Parameters:self – An instance of the Grid class
Returns:Cell size in 2D
compute_max()[source]

Compute the max values for x, y and z of the grid.

Parameters:self – An instance of the Grid class
Returns:Max values of x,y and z
cs2

Some helpful things for getting grid indices

del_cells()[source]

Delete the cells for a cells grid

Parameters:self – An instance of the Grid class
Returns:
del_lx()[source]

Delete the x-dimension of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in x-direction
del_ly()[source]

Delete the y-dimension of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in y-direction
del_lz()[source]

Delete the z-dimension of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in z-direction
del_points()[source]

Delete the points of a points grid

Parameters:self – An instance of the Grid class
Returns:
get_cells()[source]

Update the number of cells for a cells grid

Parameters:self – An instance of the Grid class
Returns:Number of cells for x, y and z
get_lx()[source]

Provide as output a tuple containing the x-size of a grid. Useful for the implementation of ‘property’.

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in x-direction
get_ly()[source]

Provide as output a tuple containing the y-size of a grid. Useful for the implementation of ‘property’.

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in y-direction
get_lz()[source]

Provide as output a tuple containing the z-size of a grid. Useful for the implementation of ‘property’.

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in z-direction
get_points()[source]

Update the number of points for a points grid

Parameters:self – An instance of the Grid class
Returns:Number of points for a points grid
idx_z(zval)[source]

Return the index of an elevation for the k-dimension of a 3D array

Parameters:zval – Elevation value
Returns:iz (int) - Index of elevation
lx

‘size’ along x of the grid.

ly

‘size’ along y of the grid.

lz

‘size’ along z of the grid.

meshup(ind='ij')[source]

Create a meshgrid representation of the grid

Parameters:self – An instance of the Grid class
Returns:A tuple containing the x,y,z-coordinates of the grid
meshup2d(ind='ij')[source]

Create a 2D meshgrid representation of the grid

Parameters:self – An instance of the Grid class
Returns:A tuple containing the x,y-coordinates of the grid
origin()[source]

To print out the origin of the grid as a tuple

Parameters:self – An instance of the Grid class
Returns:A tuple containing the origin defined in the grid
points

Number of points

print_intervals(axis='xyz')[source]

Print the intervals that constitute the simulation domain in a format like:

[ ox, ox+nx*dx] [ oy, oy+ny*dy] [ oz, oz+nz*dz]

where ox is the origin, nx is the number of points and dx is the delta between points (idem for y and z).

Parameters:axis (string) – containing [‘x’,’y’,’z’], [optional] If the default value “xyz” is used, then all the intervals are printed.
Returns:print intervals that constitute the simulation domain
set_cells(val=None)[source]

Set the number of cells for a cells grid

Parameters:self – An instance of the Grid class
Returns:Number of cells for x, y and z
set_lx(val=None)[source]

Set the x-size of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in x-direction
set_ly()[source]

Set the y-size of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in y-direction
set_lz()[source]

Set the z-size of a grid

Parameters:self – An instance of the Grid class
Returns:Tuple containing the size of a grid in z-direction
set_points(val=None)[source]

Set the number of points for a points grid

Parameters:self – An instance of the Grid class
Returns:Number of points for a points grid
shape()[source]

To print out the shape of the grid as a tuple

Parameters:
  • self – An instance of the Grid class
  • gtype (string) – ‘points’ or ‘cells’. String to decide to print the shape in terms of points or in terms of cells.
Returns:

A tuple containing the shape defined in the grid

spacing()[source]

To print out the spacing of the grid as a tuple

Parameters:self – An instance of the Grid class
Returns:A tuple containing the spacing defined in the grid
vec()[source]

Create vectors of spatial coordinates

Parameters:self – An instance of the Grid class
Returns:xv, yv, zv - Vectors of spatial coordinates
vec_node()[source]

Create vectors of spatial coordinates of bounding nodes

Parameters:self – An instance of the Grid class
Returns:xv, yv, zv - Vectors of spatial coordinates of bounding nodes
vec_x()[source]

Create vector of spatial x-coordinate

Parameters:self – An instance of the Grid class
Returns:xv - Vector of spatial x-coordinate
vec_y()[source]

Create vector of spatial y-coordinate

Parameters:self – An instance of the Grid class
Returns:yv - Vector of spatial y-coordinate
vec_z()[source]

Create vector of spatial z-coordinate

Parameters:self – An instance of the Grid class
Returns:zv - Vector of spatial z-coordinate

hyvr.sim module

Hydrogeological Virtual Reality simulation package.

Hydrogeological virtual reality (HYVR) simulator for object-based modelling of sedimentary structures

Author:Jeremy P. Bennett
Notes:Grid nodes are cell-centred!
hyvr.sim.angle(v1, v2)[source]

Return angle between two vectors in [°]

Parameters:
  • v1 – Vector 1
  • v2 – Vector 2
Returns:

angle value (float) - Angle between v1 and v2

hyvr.sim.channel_checker(param_file, ae_name, no_extpar=1, dist=0)[source]

extruded parabola checker function for quickly assessing the shape of extruded parabola inputs

Parameters:
  • param_file (str) – Parameter file location
  • ae_name (str) – Name of architectural element
  • no_extpar (float) – Number of extruded parabolas
  • dist (float) – Distance to generate extruded parabolas - defaults to mg.lx
Returns:

Plots showing shape of Ferguson extruded parabolas

hyvr.sim.curve_interp(xc, yc, spacing)[source]

Interpolate evenly spaced points along a curve. This code is based on code in an answer posted by ‘Unutbu’ on http://stackoverflow.com/questions/19117660/how-to-generate-equispaced-interpolating-values (retrieved 17/04/2017)

Parameters:
  • xc – x coordinates of curve
  • yc – y coordinates of curve
  • spacing – Spacing between points
Returns:

  • xn - x coordinates of interpolated points
  • yn - y coordinates of interpolated points

hyvr.sim.dip_rotate(azimuth_in, dip_in)[source]

Rotate dip angle based on azimuth Note that inputs and outputs are in degrees

Parameters:
  • azimuth_in – Azimuth input angle
  • dip_in – Dipping input angle
Returns:

dip_out - Azimuth output angle

hyvr.sim.dip_sets(mg, aep, znow, curve=[], select=[], azimuth_z=0)[source]

Generate dip angles and assign to the dip matrix

Parameters:
  • mg – Mesh grid object class
  • aep – Architectural element parameters (dict)
  • cruve – Tuple of x,y coordinates of extruded parabola (omitted for linear flows) - x, y coordinates of extruded parabola - vx, vy of extruded parabola flow
  • select – Model grid nodes to assign
Returns:

  • dip_out - Array of assigned dip values
  • fac_out - Array of assigned hydrofacies

hyvr.sim.ellipsoid_gradient(x, y, z, a, b, c, alpha, select, tr)[source]

Calculate dip and strike values in rotated ellipsoids

Parameters:
  • y, z (x,) – Distances to centre of ellipsoid
  • b, c (a,) – Majox/minor axes of ellipsoid
  • alpha – Rotation of ellipsoid from mean flow direction
Returns:

  • dip_g - Dipping in 3D
  • azimuth_g - Azimuth in 3D

hyvr.sim.facies(run, model, strata, hydraulics, flowtrans, elements, mg)[source]

Generate hydrofacies fields

Parameters:
  • run – Model run parameters like runname, rundir, l_dataoutputs, l_modeloutputs, etc.
  • model – Model domain parameters
  • strata – Details about the strata
  • hydraulics – Details about the hydraulics
  • flowtrans – Flow & transport simulation parameters
  • elements – Architectural elements and parameters
  • mg – Mesh grid object class
Returns:

Tuple containing:

  • probs (dict): Contains data of architectural element units and associated hydrofacies
  • params (tuple) - Contains parameters of model domain, strata, hydraulics, etc.
    • run (dict) - Model run parameters
    • model (dict) - Model domain parameters
    • strata (dict) - Strata parameters
    • hydraulics (dict) - Hydraulic properties parameters
    • flowtrans (dict) - Flow & transport simulation parameters
    • elements (dict) - Architectural elements and parameters
    • mg - Mesh grid object class
    • ae_lu - Architectural element lookup table

Return type:

(tuple)

hyvr.sim.ferguson_curve(mg, h, k, ds, eps_factor, dist=0, disp=False, ch_start=[])[source]

Simulate extruded parabola centrelines using the Ferguson (1976) disturbed meander model Implementation of AR2 autoregressive model http://onlinelibrary.wiley.com/doi/10.1002/esp.3290010403/full

Parameters:
  • mg (object class) – Mesh grid object class
  • h (float) – Height
  • k (float) – Wave number
  • ds (float) – Curve distance for calculations
  • eps_factor (float) – Random background noise
  • dist (float) – Distance to generate curves - defaults to mg.lx
  • disp (bool) – Creating display extruded parabola - extruded parabola begins at (0,0)
  • ch_start (tuple) – Starting location of channel (x,y coordinates)
Returns:

Simulated extruded parabola centerlines: storage array containing values for x coordinate, y coordinate, vx and vy

Return type:

outputs (float array)

hyvr.sim.ferguson_theta(s, eps_factor, k, h)[source]

Calculate curve direction angle

Parameters:
  • s – Steps within generated curve distance
  • eps_factor – Random background noise
  • k – Wave number
  • h – Height
Returns:

th_store (array) - Curve direction angle

hyvr.sim.gen_extpar(ch_par, mg, model, ssm, ae_array, count, ani=True)[source]
Generate extruded parabola geometries:
  • Flow regime is assumed to be reasonably constant so the major geometry of the extruded parabolas doesn’t change so much
  • ‘Migration’ of the extruded parabolas according to a shift vector
Parameters:
  • ch_par – Extruded parabola parameters
  • mg – Mesh grid object class
  • model (dict) – Model domain parameters
  • ssm (dict) – Strata parameters
  • ae_array – Array with architectural element unit details
  • count (int) – Material number and/or identifier
  • ani (bool) – Boolean if anisotropy is to be generated
  • (z_in – starting depth)
  • (thickness – Thickness of architectural element)
Returns:

Tuple containing:

  • probs (dict): Contains data of architectural element units and associated hydrofacies
    • mat - Material values
    • azim - Azimuth angles (ani == True)
    • dip - Dipping angles (ani == True)
    • fac - Facies values
    • ae_arr_i - Array with architectural element unit details
  • count (int) - Material number and/or identifier

Return type:

(tuple)

hyvr.sim.gen_sheet(sh, mg, ae_i, ae_array, count, ani=True)[source]

Generate gravel sheet with internal heterogeneity

Parameters:
  • sh – Sheet parameters
  • mg – Model grid class
  • ae_i – Architectural element lookup details [system number, z_bottom, z_top, architectural element, geometry]
  • ae_array – Architectural element array
  • count (int) – Material number and/or identifier
  • ani (bool) – Boolean if anisotropy is to be generated
Returns:

  • probs (dict) - Contains data of architectural element units and associated hydrofacies (e.g. values of azimuth, material, dipping, etc.)
  • count (int) - Material number and/or identifier

hyvr.sim.gen_trough(tr, mg, model, ae, ae_arr, count, ani=True)[source]

Create trough shapes

Parameters:
  • tr (dict) – Trough parameters
  • mg (grid class) – Model grid
  • ae (list) – Architectural element unit details
  • ae_arr (ndarray) – 3D array of system numbers
  • count (int) – Material number and/or identifier
  • ani (bool) – Boolean if anisotropy is to be generated
Returns:

  • probs (numpy array) - Grid properties
  • count (int) - Material number and/or identifier

hyvr.sim.heterogeneity(props, params)[source]

Generate internal heterogeneity

Parameters:
  • probs (list) – Data of architectural element units and associated hydrofacies (e.g. values of azimuth, material, dipping, etc.)
  • params (list) – Parameters of model domain, system, hydraulics, etc.
Returns:

  • probs (list) - Data of architectural element units and associated hydrofacies
  • params (list) - Input parameters, assigned with heterogenity

hyvr.sim.planepoint(dip_norm, x_dip, y_dip, znow, xtemp, ytemp, ztemp, select=[])[source]

Compute number of planes

Parameters:
  • dip_norm
  • x_dip – X coordinates of points on dip planes
  • y_dip – Y coordinates of points on dip planes
  • znow – Current coordinates of Z, needed to compute Z coordinates of points on dip planes
  • xtemp – X dimension of model grid nodes
  • ytemp – Y dimension of model grid nodes
  • ztemp – Z dimension of model grid nodes
  • select – Model grid nodes to consider
Returns:

set_no - Number of planes with selected model grid nodes

hyvr.sim.prob_choose(choices, probs)[source]

Get random values of an architectural element

Parameters:
  • choices – Fixed number of choices
  • probs – Contains data of architectural element units and associated hydrofacies
Returns:

choice (int) - Random value of architectural elements

hyvr.sim.rand_trough(tr, mg=False, ztr=[])[source]

Randomly generate ellipsoid geometry parameters:

Parameters:
  • tr – Ellipsoid parameters
  • mg – Meshgrid object
  • ztr – Elevation of ellipsoid centre point
Returns:

a, b, c - Length, width and depth of ellipsoid

hyvr.sim.run(param_file, flag_ow=None)[source]

Main function for HYVR generation

Parameters:
  • param_file (str) – Parameter file location
  • flag_ow – bool If True, an existing run directory is overwritten.
Returns:

Save data outputs as parameter file

hyvr.sim.save_arrays(arr_size, bg=None, mat_count=0, ani=True)[source]

Generate arrays for material properties storage

Parameters:
  • arr_size – Size of array
  • bg – List of background values for each array
  • ani (bool) – Boolean if anisotropy is to be generated
Returns:

  • ha_arr - Material values
  • fac - Facies values

hyvr.sim.save_models(realdir, realname, mg, outputs, flowtrans, props)[source]

Save HYVR outputs to standard modelling codes

Parameters:
  • run (dict) – Model run parameters
  • mg – Mesh grid object class
  • flowtrans (dict) – Flow & transport simulation parameters
  • props (dict) – Property fields (numpy arrays in dict)
Returns:

Save data outputs as .mf (MODFLOW) or .hgs (HydroGeoSphere)

hyvr.sim.save_outputs(realdir, realname, outputs, mg, outdict, suffix=None)[source]

Save data arrays to standard formats

Parameters:
  • realdir (str) – File path to save to
  • realname (str) – File name
  • outputs (str) – String codes for what type of outputs should be saved
  • mg – Mesh grid object class
  • outdict – Output directory
  • suffix (str) – Suffix for file names
Returns:

Save data outputs as .vtk (Paraview), .mat (Matlab) or .dat (Python pickle output)

hyvr.sim.scale_rotate(x, y, z, alpha=0, a=1, b=1, c=1)[source]

Scale and rotate three-dimensional trough

Parameters:
  • y, z (x,) – Spatial coordinates
  • alpha (float) – Rotation angle about the z-axis
  • b, c (a,) – Axis lengths in x, y, z directions (ellipsoid length, width, depth)
Returns:

  • select - Grid cells within ellipsoid
  • R2 - Grid of scaled and rotated values

hyvr.sim.thetaAR2(t1, t2, k, h, eps)[source]

Implementation of AR2 autoregressive model (Ferguson, 1976, Eq.15) http://onlinelibrary.wiley.com/doi/10.1002/esp.3290010403/full

Parameters:
  • t1 – theta(i-1)
  • t2 – theta(i-2)
  • k – Wavenumber
  • h – Height
  • eps – Random background noise
Returns:

2nd-order autoregression (AR2)

hyvr.utils module

Some utility functions for HFM modelling

Authors:Jeremy P. Bennett, with help from Alessandro Comunian and Samuel Scherrer
Notes:
hyvr.utils.calc_norm(x)[source]

Calculate norm (compute the complex conjugate from ‘x’)

Parameters:x – Input parameter
Returns:Complex conjugate of x
hyvr.utils.dem_load(fn)[source]

Load data from ESRI-style ASCII-file.

Parameters:fn (str) – Directory and file name for save
Returns:
  • data (numpy array) - Data from ERSI-style ASCII-file
  • meta (dict) - Dict with grid metadata
hyvr.utils.dem_save(fn, data, gro)[source]

Save DEM data to ESRI-style ASCII-file

Parameters:
  • fn (str) – Directory and file name for save
  • data (numpy array) – DEM data
  • gr (object class) – grid.Grid() object class
Returns:

Save DEM data to ESRI-style ASCII-file

hyvr.utils.load_gslib(fn)[source]

Load .gslib files. This has been appropriated from the HPGL library https://github.com/hpgl/hpgl/blob/master/src/geo_bsd/routines.py commit b980e15ad9b1f7107fd4fa56ab117f45553be3aa

Parameters:fn (str) – .gslib file path and name
Returns:gslib_dict (dict) - properties
hyvr.utils.load_pickle(pickfile)[source]

Pickle input file

Parameters:pickfile – Input file
Returns:data (dict) - Pickled data of input file
hyvr.utils.matlab_save(fn, data)[source]

Save numpy array to .mat file for use in matlab.

Parameters:
  • fn (str) – File name (ending with .mat)
  • data (numpy array) – Data to save
Returns:

Save a dictionary of names and arrays into a MATLAB-style .mat file. This saves the array objects in the given dictionary to a MATLAB- style .mat file.

hyvr.utils.mf6_vtr(fhead, mg, fout)[source]

Convert a MODFLOW 6 Binary head file into vtr suitable for visualisation in ParaView

hyvr.utils.read_lu(sq_fp)[source]

Load user-defined strata (architectural element lookup table), split the data based on a delimiter and return it as a new list

Parameters:sq_fp – Load user-defined strata (architectural element lookup table)
Returns:-Values of architectural element lookup table
Return type:ssm_lu (list)
hyvr.utils.rotate_ktensor(count, aniso, azimuth, dip, k_in)[source]

Create a rotated K tensor

Parameters:
  • count (int) – Material number and/or identifier
  • aniso – Anisotropy
  • azimuth – Azimuth angles
  • dip – Dipping angles
  • k_in
Returns:

k_rotate - Rotated K tensor

hyvr.utils.round_x(x, base=1, prec=2)[source]

Round to the nearest z-increment (Refer to http://stackoverflow.com/questions/2272149/round-to-5-or-other-number-in-python)

Parameters:
  • x (float) – Input parameter
  • base (int) – Base parameter for avoiding floating-point values
  • prec – Precision of rounding
Returns:

Rounded value of nearest z-increment

hyvr.utils.specsim(gr, var, corl, twod=False, covmod='gau')[source]

Generate random variables stationary covariance function using spectral techniques of Dietrich & Newsam (1993)

Parameters:
  • gr – Grid class object
  • var – Variance
  • corl – Tuple of correlation length of random variable
  • twod – Flag for two-dimensional simulation
  • covmod – Which covariance model to use (‘gau’ = Gaussian, ‘exp’ = Exponential).
Returns:

bigy - Random gaussian variable. Real part of a complex array, created via inverse DFT

hyvr.utils.to_hgs(hgspath, mg, flowtrans, ktensors, poros)[source]

Convert HYVR outputs to HydroGeoSphere inputs

Parameters:
  • hgspath (str) – Path where to save HGS output file
  • ktensors – Array with tensor values of K
  • poros – Array with values of porosity
Returns:

  • val_fmts (dict) - Dictionary with values of K tensor and porosity
  • val_filepath - file name of HGS output file

hyvr.utils.to_mf6(mfdir, runname, mg, flowtrans, props)[source]

Convert HYVR outputs to MODFLOW6 inputs

Parameters:
  • mfdir – Directory of MODFLOW model object
  • runname – Run name
  • mg – Mesh grid object class
  • flowtrans (dict) – Flow & transport simulation parameters
  • props (dict) – Parameter fields as numpy arrays (k_iso: Hydraulic conductivity of HYVR, anirat: Background anistropic ratio (K_h/K_v anisotropy ratio))
Returns:

  • mf - MODFLOW model object
  • dis - Discretization of modflow object
  • bas - BAS package of modflow model
  • lpf - LPF package of modflow model
  • oc - OC package of modflow model
  • pcg - pcg package of modflow model

hyvr.utils.to_modflow(mfdir, mg, flowtrans, props)[source]

Convert HYVR outputs to MODFLOW inputs

Parameters:
  • mfdir – Directory of MODFLOW model object
  • mg – Mesh grid object class
  • flowtrans (dict) – Flow & transport simulation parameters
  • props (dict) – Parameter fields as numpy arrays (anirat: Background anistropic ratio (K_h/K_v anisotropy ratio)
Returns:

  • mf - MODFLOW model object
  • dis - Discretization of modflow object
  • bas - BAS package of modflow model
  • lpf - LPF package of modflow model
  • oc - OC package of modflow model
  • pcg - pcg package of modflow model

hyvr.utils.to_vtr(data_dict, file_name, grid, points=None)[source]

Save a numpy array into a .vtr rectilinear grid of voxels using pyevtk

Parameters:
  • data_dict (numpy array) – e.g. {‘fac’: fac, ‘mat’: mat}
  • file_name (str) – Name of the file for the output.
  • grid (class Grid) – The information about the grid can be also provided as a grid object.
  • points (bool) – Set as True for exporting cell-centered points for contouring (e.g., hydraulic heads)
Returns:

A VTK STRUCTURED_POINTS dataset file containing the input numpy data.

hyvr.utils.try_makefolder(makedir)[source]

Create modflow output folder

hyvr.utils.virtual_boreholes(data_dict, d, l, file_out=None, vals=[], opts=[])[source]

Perform ‘virtual’ borehole sampling of parameter field

Parameters:
  • data_dict (dict) – Data to sample
  • d (list) – 3-tuple of model grid cell dimensions
  • l (list) – 3-tuple of total model dimensions/lengths
  • (list (file_out) – basefile path, list of file types): Output filename and path
  • vals (list) – Parameter fields to include
  • opts (dict) – Sampling options opts[‘noBH’] (int): Random sampling opts[‘grid_spacing’]: float, or list of floats [x spacing, y spacing] Grid sample spacing opts[‘grid_n’]: int, or list of ints [n in x, n in y] Number of grid nodes per x,y dimensions opt[‘lnK’] (bool): Natural logarithm of isotropic hydraulic conductivity
Returns:

Pandas DataFrame class

Return type:

bh_df

hyvr.parameters module

This is the main module for parsing the parameter file. There are two other modules strongly linked ot this one:

  • options: contains definitions of the possible options of a parameter file
  • option_parsing: provides the classes Section and Option for simpler parsing and validation.
Authors:Jeremy P. Bennett, Samuel Scherrer
hyvr.parameters.model_setup(pf)[source]

Set up model using grid.Grid() class and assign parameters

Parameters:pf (str) – Parameter file path
Returns:
  • run (dict) - Model run parameters
  • mod (dict) - Model domain parameters
  • sequences (dict) - Sequence parameters
  • hydraulics (dict) - Hydraulic properties parameters
  • flowtrans (dict) - Flow & transport simulation parameters
  • elements (dict) - Architectural elements and parameters
  • model_grid (object class) - Grid object class
hyvr.parameters.parameters(inifile)[source]

Get parameters for hierarchical facies modelling from a .ini-file

Parameters:inifile (str) – Parameter file path
Returns:
  • run (dict) - Model run parameters
  • model (dict) - Model domain parameters
  • sequences (dict) - Sequence parameters
  • hydraulics (dict) - Hydraulic properties parameters
  • flowtrans (dict) - Flow & transport simulation parameters
  • elements (dict) - Architectural elements and parameters
hyvr.parameters.set_up_directories(run, configfile, flag_ow=None)[source]

This functions creates the necessary directories (modeldir, rundir). It also stores the used ini-file in the rundir.

Parameters:
  • run (dict) – parsed run-section of the config file
  • configfile (str) – path to config file
  • flag_ow (bool) – Whether to overwrite the old run directory. If it is None, the option from the config file will be chosen instead (default in configfile: False)

hyvr.options module

This module contains all options for the ini-files. If you want to add an option, add it here and potentially do some post-processing in parameters.

Authors:Samuel Scherrer, Jeremy P. Bennett

hyvr.option_parsing module

This module contains the classes Section and Option to simplify parsing and validating ini-files.

Author:Samuel Scherrer
class hyvr.option_parsing.Section(name, options)[source]

Bases: object

parse(section_dict)[source]

Parses and validates options of the section given a dictionary containing key-value pairs of the options. If options are not present, default values might be set.

Parameters:section_dict (dictionary) –

Dictionary of section values. This can for example be obtained using configparser:

p = configparser.ConfigParser()
p.read(filename)
section_dict = dict(p[section_name])
Returns:section_dict – The same dictionary with parsed and validated values
Return type:dictionary
class hyvr.option_parsing.Option(name, dtype, optional=False, default=None, shape=-1, datatype=None, validation_func=<function Option.<lambda>>, alternatives=[])[source]

Bases: object

Parameters:
  • name (string) – name of the option
  • dtype (type) – type of the value of the option, e.g. float, str, int, or list. If dtype is list, every entry of the list must have the same type.
  • optional (bool, optional (default: False)) – whether the option is optional or not
  • default (optional, (default: None)) – if optional=True, this default value will be used if this option is not given.
  • shape (int or string or list/tuple of ints and/or strings, optional (only) –

    required for lists, default: None) If dtype is list, this is the shape of the list. There are several possibilities how to use this option:

    • if shape=n where n is a nonnegative integer, the value must be a list with length n.
    • if shape=-1 the value can have an arbitrary shape.
    • if shape="option1", the value of this option must have the same shape as “option1”. This is especially useful if the shape of “option1” is set to -1.
    • if shape=[2, 3], the value must be a list of lists, where the outermost list has length 2 and the inner lists all have length 3. This also works for more than 2 dimensions.
    • if shape=[2, -1, 3], the value must be a list of lists of lists. The outermost list must again have length 2, the innermost lists must have length 3, and the lists at the intermediate level can have any length (even different lengths).
    • if shape=[2, "option1", 3], the values must again be a list of lists similar to above, but now the lists at the intermediate level must have the same length as “option1”.
    • if shape=[2, [1, 2], [[3], [3, 3]]], the value must be a list of lists of lists. The outermost list must again have length 2. The two lists it contains have length 1 and length 2. The innermost lists all must have length 3.

    It’s also possible to only give the innermost value(s), e.g. for a list with shape=[2, 3] only the value 18. This will then be expanded to [[18, 18, 18], [18, 18, 18]]. Similarly, [18, 19, 20] would be expanded to [[18, 19, 20], [18, 19, 20]]. This expansion obviously doesn’t work if shape=-1. If shape=[2, -1, 3], expansion is possible if the given value is e.g. [[1, 2, 3], [1, 2, 3], [1, 2, 3]], but not if only [1, 2, 3] is given, since the length of the second dimension must be determined from the given value.

  • datatype (int, float or string, optional (only required for lists, default: None)) – Type of the innermost elements of the option in case the option is a list.
  • validation_func (function of one argument that returns a boolean, optional.) – Returns true if the value (for lists or lists of lists this applies to all values) is valid.
  • alternatives (list of strings) – List of alternative names for this option. If the option is not found in the given dictionary while parsing, these alternatives are used if they exist. If an alternative is found (searching happens in the supplied order), the option is stored under it’s standard name.
parse(section)[source]

Parses the option based on it’s attributes.

exception hyvr.option_parsing.MissingSectionError(sectionname)[source]

Bases: Exception

exception hyvr.option_parsing.MissingOptionError(optionname, sectionname)[source]

Bases: Exception

exception hyvr.option_parsing.ShapeError(optionname, sectionname)[source]

Bases: Exception

hyvr.option_parsing.assert_exists(option, dictionary, section_name)[source]