DPFEHM.addboundaryconditionsMethod

addboundaryconditions(hfree, dirichletnodes, dirichleths, isfreenode, nodei2freenodei)

Returns an array giving the pressures at all nodes, including the nodes with Dirichlet boundary conditions.

Arguments

  • hfree: an array giving values of h at each cell that is not constrained by a Dirichlet boundary condition
  • dirichletnodes: the indices of the cells that are constrained by a Dirichlet boundary condition
  • dirichleth: the values of h at the cells with Dirichlet boundary conditions
  • isfreenode: a boolean array such that isfreenode[i] is true if cell i does not have a Dirichlet boundary condition
  • nodei2freenodei: a dictionary that maps a node index to a free node index
DPFEHM.amg_solverMethod

amg_solver(A, b; kwargs...)

Return the solution to Ax=b using an algebraic multigrid solver, which the kwargs being passed to the AlgebraicMultigrid solver

DPFEHM.cholesky_solverMethod

cholesky_solver(A, b; kwargs...)

Return the solution to Ax=b using a Cholesky factorization solver

DPFEHM.darcy_velocityMethod

darcy_velocity(h, Ks, mins, maxs, ns)

Return a function that computes the Darcy velocity on a regular grid

Arguments

  • h: pressure
  • Ks: permeability
  • mins: the minimum x[/y/z] coordinates of the cube
  • maxs: the maximum x[/y/z] coordinates of the cube
  • ns: the number of grid points in the x, [y, and z] directions
DPFEHM.effective_saturationMethod

effective_saturation(alpha::Number, psi::Number, N::Number)

Return the effective saturation

Arguments

  • alpha: van Genucthen's α parameter
  • psi: pressure head
  • N: van Genuchten's N parameter
DPFEHM.getfreenodesMethod

getfreenodes(n, dirichletnodes)

Returns an array indicating whether or not a node is free, a dictionary that maps a global index into a free node index, and an array that maps a free node index into a global index. A free node is one that is not defined by a Dirichlet boundary coundition. That is, a free node, is one that is not listed in dirichletnodes.

Arguments

  • n: the number of cells on the mesh
  • dirichletnodes: the indices of the cells that are constrained by a Dirichlet boundary condition
DPFEHM.getpointsMethod

getpoints(min, max, n)

Return a range from min to max with n points

DPFEHM.getuItersMethod

getuIters(cval, f, nz, nx, nt, dz, dx, dt)

Return the wave amplitude at the next nt time steps starting from a zero initial condition

Arguments

  • cval: velocity field
  • f: forcing
  • nz: number of grid cells in the z direction
  • nx: number of grid cells in the x direction
  • nt: number of time steps
  • dz: distance between grid cells in the z direction
  • dx: distance between grid cells in the x direction
  • dt: time step
DPFEHM.groundwater_residualsFunction

groundwater_residuals(h, Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, Qs, specificstorage, volumes; linear_solver::Function=amg_solver, kwargs...)

Return the residuals at the free nodes of a finite volume discretization of the groundwater flow equation.

Arguments

  • h: pressure
  • Ks: permeability
  • neighbors: array of pairs indicating which cells share an interface
  • areasoverlengths: array with the same length as neighbors that gives the interfacial area divided by the length between the two cell centers
  • dirichletnodes: array of indices indicating which nodes have Dirichlet boundary conditions
  • dirichleths: array of pressures at the Dirichlet boundary (length is equal to the number of cells on the grid)
  • Qs: array of fluxes (length is equal to the number of cells on the grid)
  • specificstorage: array of the specific storage associated with each cell (length is equal to the number of cells on the grid)
  • volumes: array of the volume of each each cell (length is equal to the number of cells on the grid)
DPFEHM.groundwater_steadystateMethod

groundwater_steadystate(Ks, neighbors, areasoverlengths, dirichletnodes, dirichleths, Qs; linear_solver::Function=amg_solver, kwargs...)

Return the solution to a steady state groundwater flow problem

Arguments

  • Ks: permeability
  • neighbors: array of pairs indicating which cells share an interface
  • areasoverlengths: array with the same length as neighbors that gives the interfacial area divided by the length between the two cell centers
  • dirichletnodes: array of indices indicating which nodes have Dirichlet boundary conditions
  • dirichleths: array of pressures at the Dirichlet boundary (length is equal to the number of cells on the grid)
  • Qs: array of fluxes (length is equal to the number of cells on the grid)
DPFEHM.hmMethod

hm(x)

Returns the harmonic mean of x[1] and x[2]

DPFEHM.hmMethod

hm(x, y)

Returns the harmonic mean of x and y

DPFEHM.krMethod

kr(x)

Returns kr(x[1], x[2], x[3])

DPFEHM.krMethod

kr(psi, alpha, N)

Return the van Genuchten relative permeability

Arguments

  • psi: pressure head
  • alpha: van Genucthen's α parameter
  • N: van Genuchten's N parameter
DPFEHM.load_ugeMethod

load_uge(filename)

Return mesh information (cell coordinates, cell volumes, cell neighbors, interfacial areas, and distances between neighboring nodes) from the ASCII .uge filename

For more information on the format, see: https://documentation.pflotran.org/userguide/cards/subsurface/grids/unstructuredexplicit_grid.html

DPFEHM.make_wave2diter_pullbackMethod

make_wave2diter_pullback(args...)

Return a pullback function for wave2diter_residuals. The arguments are the same as for wave2diter_residuals

DPFEHM.regulargrid1dMethod

regulargrid1d(mins, maxs, ns, dy, dz)

Return mesh information (cell coordinates, cell volumes, cell neighbors, interfacial areas divided by the length between cell centers)

Arguments

  • mins: the minimum x coordinates of the cube
  • maxs: the maximum x coordinates of the cube
  • ns: the number of grid points in the x directions
  • dy: the size of the cells in the y dimension
  • dz: the size of the cells in the z dimension
DPFEHM.regulargrid2dMethod

regulargrid2d(mins, maxs, ns, dz)

Return mesh information (cell coordinates, cell volumes, cell neighbors, interfacial areas divided by the length between cell centers)

Arguments

  • mins: the minimum x and y coordinates of the cube
  • maxs: the maximum x and y coordinates of the cube
  • ns: the number of grid points in the x and y directions
  • dz: the size of the cells in the z dimension
DPFEHM.regulargrid3dMethod

regulargrid3d(mins, maxs, ns)

Return mesh information (cell coordinates, cell volumes, cell neighbors, interfacial areas divided by the length between cell centers)

Arguments

  • mins: the minimum x, y, and z coordinates of the cube
  • maxs: the maximum x, y, and z coordinates of the cube
  • ns: the number of grid points in the x, y, and z directions
DPFEHM.richards_residualsFunction

wave2diter(c, f, u, uim1, uim2, nz, nx, nt, dz, dx, dt)

Return the residuals at the free nodes of a 2D finite difference discretization of the wave equation

Arguments

  • c: velocity field
  • f: forcing
  • u: wave
  • uim1: wave at the previous time step
  • uim2: wave two time steps ago
  • nz: number of grid cells in the z direction
  • nx: number of grid cells in the x direction
  • nt: number of time steps
  • dz: distance between grid cells in the z direction
  • dx: distance between grid cells in the x direction
  • dt: time step
DPFEHM.richards_steadystateMethod

richards_steadystate(Ks, neighbors, areasoverlengths, dirichletnodes, dirichletpsis, coords, alphas, Ns, Qs; callback=soln->nothing, kwargs...)

Return the solution to a steady state unsaturated groundwater problem

Arguments

  • psi0: initial guess for the pressure
  • Ks: permeability
  • neighbors: array of pairs indicating which cells share an interface
  • areasoverlengths: array with the same length as neighbors that gives the interfacial area divided by the length between the two cell centers
  • dirichletnodes: array of indices indicating which nodes have Dirichlet boundary conditions
  • dirichletpsis: array of pressures at the Dirichlet boundary (length is equal to the number of cells on the grid)
  • coords: the coordinates of the cell centers
  • alphas: van Genuchten α parameters for each cell
  • Ns: van Genuchten N parameters for each cell
  • Qs: array of fluxes (length is equal to the number of cells on the grid)
DPFEHM.transport_residualsFunction

transport_residuals(u, vxs, vys, vzs, Ds, neighbors, areasoverlengths, dirichletnodes, dirichletus, Qs, volumes, coords; use_upwind=true)

Return the residuals of a finite volume discretization of the advection-dispersion equation.

#Arguments

  • u: the concentration of the species undergoing transport
  • vxs: the x components of the velocity at each cell
  • vys: the y components of the velocity at each cell
  • vzs: the z components of the velocity at each cell
  • Ds: the dispersion coefficient at each cell
  • neighbors: array of pairs indicating which cells share an interface
  • areasoverlengths: array with the same length as neighbors that gives the interfacial area divided by the length between the two cell centers
  • dirichletnodes: array of indices indicating which nodes have Dirichlet boundary conditions
  • dirichletus: array of concentrations at the Dirichlet boundary (length is equal to the number of cells on the grid)
  • Qs: array of fluxes (length is equal to the number of cells on the grid)
  • volumes: array of the volume of each each cell (length is equal to the number of cells on the grid)
DPFEHM.uoneiterMethod

uoneiter(ui, uim1, uim2, c, f, nz, nx, nt, dz, dx, dt)

Return the wave amplitude at the next time step

Arguments

  • ui: wave at the current time step
  • uim1: wave at the previous time step
  • uim2: wave two time steps ago
  • c: velocity field
  • f: forcing
  • nz: number of grid cells in the z direction
  • nx: number of grid cells in the x direction
  • dz: distance between grid cells in the z direction
  • dx: distance between grid cells in the x direction
  • dt: time step