API
Initialization and configuration
HYPRE.Init
— FunctionInit(; 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.
Matrix/vector creation
HYPRE.start_assemble!
— FunctionHYPRE.start_assemble!(A::HYPREMatrix) -> HYPREMatrixAssembler
+API · HYPRE.jl API
Initialization and configuration
HYPRE.Init
— FunctionInit(; 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.
sourceMatrix/vector creation
HYPRE.start_assemble!
— FunctionHYPRE.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!
.
sourceHYPRE.assemble!
— FunctionHYPRE.assemble!(A::HYPREMatrixAssembler, i, j, a::Matrix)
+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!
.
sourceHYPRE.assemble!
— FunctionHYPRE.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!
.
sourceHYPRE.finish_assemble!
— FunctionHYPRE.finish_assemble!(A::HYPREMatrixAssembler)
+A.b[i] += b
See also: HYPRE.start_assemble!
, HYPRE.finish_assemble!
.
sourceHYPRE.finish_assemble!
— FunctionHYPRE.finish_assemble!(A::HYPREMatrixAssembler)
HYPRE.finish_assemble!(A::HYPREVectorAssembler)
-HYPRE.finish_assemble!(A::HYPREAssembler)
Finish the assembly. This synchronizes the data between processors.
sourceSolvers and preconditioners
HYPRE.solve!
— Functionsolve!(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
.
sourceHYPRE.solve
— Functionsolve(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!
.
sourceHYPRE.HYPRESolver
— TypeHYPRESolver
Abstract super type of all the wrapped HYPRE solvers.
sourceHYPRE.BiCGSTAB
— TypeBiCGSTAB(; settings...)
Create a BiCGSTAB
solver. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.BoomerAMG
— TypeBoomerAMG(; settings...)
Create a BoomerAMG
solver/preconditioner. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.FlexGMRES
— TypeFlexGMRES(; settings...)
Create a FlexGMRES
solver. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.GMRES
— TypeGMRES(; settings...)
Create a GMRES
solver. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.Hybrid
— TypeHybrid(; settings...)
Create a Hybrid
solver. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.ILU
— TypeILU(; settings...)
Create a ILU
solver/preconditioner. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.PCG
— TypePCG(; settings...)
Create a PCG
solver. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.ParaSails
— TypeParaSails(comm=MPI.COMM_WORLD; settings...)
Create a ParaSails
preconditioner. See HYPRE API reference for details and supported settings.
External links
sourceHYPRE.GetNumIterations
— FunctionHYPRE.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
.
sourceHYPRE.GetFinalRelativeResidualNorm
— FunctionHYPRE.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
.
sourceSettings
This document was generated with Documenter.jl version 1.7.0 on Sunday 29 September 2024. Using Julia version 1.10.5.
+HYPRE.finish_assemble!(A::HYPREAssembler)
Finish the assembly. This synchronizes the data between processors.
Solvers and preconditioners
HYPRE.solve!
— Functionsolve!(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
.
HYPRE.solve
— Functionsolve(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!
.
HYPRE.HYPRESolver
— TypeHYPRESolver
Abstract super type of all the wrapped HYPRE solvers.
HYPRE.BiCGSTAB
— TypeBiCGSTAB(; settings...)
Create a BiCGSTAB
solver. See HYPRE API reference for details and supported settings.
External links
HYPRE.BoomerAMG
— TypeBoomerAMG(; settings...)
Create a BoomerAMG
solver/preconditioner. See HYPRE API reference for details and supported settings.
External links
HYPRE.FlexGMRES
— TypeFlexGMRES(; settings...)
Create a FlexGMRES
solver. See HYPRE API reference for details and supported settings.
External links
HYPRE.GMRES
— TypeGMRES(; settings...)
Create a GMRES
solver. See HYPRE API reference for details and supported settings.
External links
HYPRE.Hybrid
— TypeHybrid(; settings...)
Create a Hybrid
solver. See HYPRE API reference for details and supported settings.
External links
HYPRE.ILU
— TypeILU(; settings...)
Create a ILU
solver/preconditioner. See HYPRE API reference for details and supported settings.
External links
HYPRE.PCG
— TypePCG(; settings...)
Create a PCG
solver. See HYPRE API reference for details and supported settings.
External links
HYPRE.ParaSails
— TypeParaSails(comm=MPI.COMM_WORLD; settings...)
Create a ParaSails
preconditioner. See HYPRE API reference for details and supported settings.
External links
HYPRE.GetNumIterations
— FunctionHYPRE.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
.
HYPRE.GetFinalRelativeResidualNorm
— FunctionHYPRE.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
.