Class FAMesh#
Public functions#
- FAMesh.flowAccumulation()[source]#
This function is the main entry point for flow accumulation computation.
Note
Flow accumulation (FA) calculations are a core component of landscape evolution models as they are often used as proxy to estimate flow discharge, sediment load, river width, bedrock erosion, and sediment deposition. Until recently,
goSPL model computes the flow discharge from FA and the net precipitation rate using a parallel implicit drainage area (IDA) method proposed by Richardson et al., 2014 but adapted to unstructured grids.
It calls the following private functions:
_buildFlowDirection
_solve_KSP
_distributeDownstream
- FAMesh.matrixFlow(flowdir, dep=None)[source]#
This function defines the flow direction matrices.
Note
The matrix is built incrementally looping through the number of flow direction paths defined by the user. It proceeds by assembling a local Compressed Sparse Row (CSR) matrix to a global PETSc matrix.
When setting up the flow matrix in PETSc, we preallocate the non-zero entries of the matrix before starting filling in the values. Using PETSc sparse matrix storage scheme has the advantage that matrix-vector multiplication is extremely fast.
The matrix coefficients consist of weights (comprised between 0 and 1) and calculated based on the number of downslope neighbours and proportional to the slope.
- Parameters:
flow – boolean to compute matrix for either downstream water or sediment transport
dep – deposition flux coefficient in case where the sediment transport/deposition term is considered.
Private functions#
- FAMesh._buildFlowDirection(h, down=True)[source]#
This function builds from neighbouring slopes the flow directions. It calls a fortran subroutine that locally computes for each vertice:
the indices of receivers (downstream) nodes depending on the desired number of flow directions (SFD to MFD).
the distances to the receivers based on mesh resolution.
the associated weights calculated based on the number of receivers and proportional to the slope.
- Parameters:
h – elevation in the form of a Numpy Array
down – boolean to indicate whether the filled elevation needs to be considered or not.
- FAMesh._distributeDownstream(pitVol, FA, hl, step, ice=False)[source]#
In cases where rivers flow in depressions, they might fill the sink completely and overspill or remain within the depression, forming a lake. This function computes the excess of water (if any) able to flow dowstream.
Important
The excess water is then added to the downstream flow accumulation (FA) and used to estimate rivers’ erosion.
- Parameters:
pitVol – volume of depressions
FA – excess flow accumulation array
hl – current elevation array
step – downstream distribution step
ice – boolean indicating where the ice flow is considered or not.
- Returns:
pitVol, excess, nFA (updated volume in each depression, boolean set to True is excess flow remains to be distributed and new flow accumulation values)
- FAMesh._matrix_build(nnz=(1, 1))[source]#
Creates a sparse PETSc matrix.
Note
To achieve good performance during matrix assembly, the function preallocates the matrix storage.
- Parameters:
nnz – array containing the number of nonzeros in the various rows
- Returns:
sparse PETSc matrix
- FAMesh._matrix_build_diag(V, nnz=(1, 1))[source]#
Builds a PETSc diagonal matrix based on a given array V
- Parameters:
V – diagonal data array
nnz – array containing the number of nonzero blocks
- Returns:
sparse PETSc matrix
- FAMesh._solve_KSP(guess, matrix, vector1, vector2)[source]#
PETSc scalable linear equations solvers (KSP) component provides Krylov subspace iterative method and a preconditioner. Here, flow accumulation solution is obtained using PETSc Richardson solver (richardson) with block Jacobian preconditioning (bjacobi).
Note
The solver choice was made based on the convergence results from Richardson et al. (2014) but can be changed if better solver and preconditioner combinations are found.
Using such iterative method allows for an initial guess to be provided. When this initial guess is close to the solution, the number of iterations required for convergence dramatically decreases. Here the flow discharge solution from previous time step can be passed as an initial guess to the solver as discharge often exhibits little change between successive time intervals.
- Parameters:
guess – Boolean specifying if the iterative KSP solver initial guess is nonzero (when provided it corresponds to the previous flow discharge values).
matrix – PETSc sparse matrix used by the KSP solver composed of diagonal terms set to unity (identity matrix) and off-diagonal terms (weights between 0 and 1). The weights are calculated based on the number of downslope neighbours (based on the chosen number of flow direction directions) and are proportional to the slope.
vector1 – PETSc vector corresponding to the local volume of water available for runoff during a given time step (e.g. voronoi area times local precipitation rate).
vector2 – PETSc vector corresponding to the unknown flow discharge values.
- Returns:
vector2 PETSc vector of the new flow discharge values