Class FAMesh#

class flow.flowplex.FAMesh(*args, **kwargs)[source]#

This class calculates drainage area in an implicit, iterative manner using PETSc solvers. It accounts for multiple flow direction paths (SFD to MFD) based on user input declaration.

Note

The class follows the parallel approach described in Richardson et al., 2014.

Methods

flowAccumulation()

This function is the main entry point for flow accumulation computation.

matrixFlow(flowdir[, dep])

This function defines the flow direction matrices.

Initialise

__init__(*args, **kwargs)

The initialisation of FAMesh class consists in the declaration of PETSc vectors and matrices.

Public Methods

matrixFlow(flowdir[, dep])

This function defines the flow direction matrices.

flowAccumulation()

This function is the main entry point for flow accumulation computation.

Private Methods

_matrix_build([nnz])

Creates a sparse PETSc matrix.

_matrix_build_diag(V[, nnz])

Builds a PETSc diagonal matrix based on a given array V

_make_reasons(reasons)

Provides reasons for PETSc error...

_solve_KSP2(matrix, vector1, vector2)

Solution of Krylov subspace iterative method (PETSc scalable linear equations solvers - KSP) implemented using the Flexible Generalized Minimal Residual method (fgmres) with Additive Schwarz preconditioning (asm).

_solve_KSP(guess, matrix, vector1, vector2)

PETSc scalable linear equations solvers (KSP) component provides Krylov subspace iterative method and a preconditioner.

_buildFlowDirection(h[, down])

This function builds from neighbouring slopes the flow directions.

_potentialLakeEvap()

Per-pit max evaporation volume (m^3) for one timestep, assuming the lake fills to the spillover (= max-fill surface).

_distributeDownstream(pitVol, FA, hl, step)

In cases where rivers flow in depressions, they might fill the sink completely and overspill or remain within the depression, forming a lake.

Public functions#

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.

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:

  1. _buildFlowDirection

  2. _solve_KSP

  3. _distributeDownstream

Private functions#

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._make_reasons(reasons)[source]#

Provides reasons for PETSc error…

FAMesh._solve_KSP2(matrix, vector1, vector2)[source]#

Solution of Krylov subspace iterative method (PETSc scalable linear equations solvers - KSP) implemented using the Flexible Generalized Minimal Residual method (fgmres) with Additive Schwarz preconditioning (asm).

Note

This function is used if the KSP convergence failed using the PETSc Richardson solver with block Jacobian preconditioning.

Parameters:
  • 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 (i.e., user-defined number of flow 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

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

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._potentialLakeEvap()[source]#

Per-pit max evaporation volume (m^3) for one timestep, assuming the lake fills to the spillover (= max-fill surface).

Surface = Σ larea over cells with pitIDs >= 0 (the lake’s potential extent at spillover). The inIDs == 1 mask prevents MPI halo double-counting; the per-pit total is then Allreduce-summed across ranks so every rank ends up with the same evap budget array.

See DESIGN_EVAPORATION.md §1 D3 for the design rationale (we use max-fill rather than current-fill to avoid the circular dependency between fill level and evap rate).

Returns:

numpy array of shape (nbpits,) with per-pit evap m^3

FAMesh._distributeDownstream(pitVol, FA, hl, step)[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

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)