API

Initialization and configuration

HYPRE.InitFunction
Init(; finalize_atexit=true)

Wrapper around HYPRE_Init. If finalize_atexit is true a Julia exit hook is added, which calls HYPRE_Finalize. This method will also call MPI.Init unless MPI is already initialized.

Note: This function must be called before using HYPRE functions.

source

Matrix/vector creation

HYPRE.start_assemble!Function
HYPRE.start_assemble!(A::HYPREMatrix)                 -> HYPREMatrixAssembler
HYPRE.start_assemble!(b::HYPREVector)                 -> HYPREVectorAssembler
HYPRE.start_assemble!(A::HYPREMatrix, b::HYPREVector) -> HYPREAssembler

Initialize a new assembly for matrix A, vector b, or for both. This zeroes out any previous data in the arrays. Return a HYPREAssembler with allocated data buffers needed to perform the assembly efficiently.

See also: HYPRE.assemble!, HYPRE.finish_assemble!.

source
HYPRE.assemble!Function
HYPRE.assemble!(A::HYPREMatrixAssembler, i, j, a::Matrix)
HYPRE.assemble!(A::HYPREVectorAssembler, i,    b::Vector)
HYPRE.assemble!(A::HYPREAssembler,       ij,   a::Matrix, b::Vector)

Assemble (by adding) matrix contribution a, vector contribution b, into the underlying array(s) of the assembler at global row indices i and column indices j.

This is roughly equivalent to:

# A.A::HYPREMatrix
A.A[i, j] += a

# A.b::HYPREVector
A.b[i] += b

See also: HYPRE.start_assemble!, HYPRE.finish_assemble!.

source
HYPRE.finish_assemble!Function
HYPRE.finish_assemble!(A::HYPREMatrixAssembler)
HYPRE.finish_assemble!(A::HYPREVectorAssembler)
HYPRE.finish_assemble!(A::HYPREAssembler)

Finish the assembly. This synchronizes the data between processors.

source

Solvers and preconditioners

HYPRE.solve!Function
solve!(solver::HYPRESolver, x::HYPREVector, A::HYPREMatrix, b::HYPREVector)

Solve the linear system A x = b using solver with x as the initial guess. The approximate solution is stored in x.

See also solve.

source
HYPRE.solveFunction
solve(solver::HYPRESolver, A::HYPREMatrix, b::HYPREVector) -> HYPREVector

Solve the linear system A x = b using solver and return the approximate solution.

This method allocates an initial guess/output vector x, initialized to 0.

See also solve!.

source
HYPRE.GetNumIterationsFunction
HYPRE.GetNumIterations(s::HYPRESolver)

Return number of iterations during the last solve with solver s.

This function dispatches on the solver to the corresponding C API wrapper LibHYPRE.HYPRE_$(Solver)GetNumIterations.

source
HYPRE.GetFinalRelativeResidualNormFunction
HYPRE.GetFinalRelativeResidualNorm(s::HYPRESolver)

Return the final relative residual norm from the last solve with solver s.

This function dispatches on the solver to the corresponding C API wrapper LibHYPRE.HYPRE_$(Solver)GetFinalRelativeResidualNorm.

source