Module giving access to the HYPRE solvers for sparse linear systems.
Note 1. Solver naming convention: X_Y_Z, where X: HYPRE interface (S: Struct, SS: SStruct, IJ: IJ), Y: Solver (e.g. PCG), Z: Preconditioner (e.g. DS). Available solvers: S_SMG_NONE, S_PFMG_NONE, S_PCG_NONE, IJ_AMG_NONE, IJ_PCG_DS, IJ_PCG_AMG
Note 2. GPU support is enabled with IJ_AMG_NONE. To be able to use GPU acceleration, LEAP must be compiled against GPU-enabled HYPRE.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | private, | parameter | :: | HYPRE_INT_IJ | = | 1 | |
integer, | private, | parameter | :: | HYPRE_INT_SStruct | = | 3 | |
integer, | private, | parameter | :: | HYPRE_INT_Struct | = | 2 | |
character(len=*), | private, | parameter | :: | HYPRE_SOL_IJ_AMG_NONE | = | 'IJ-AMG-NONE' | |
character(len=*), | private, | parameter | :: | HYPRE_SOL_IJ_PCG_AMG | = | 'IJ-PCG-AMG' | |
character(len=*), | private, | parameter | :: | HYPRE_SOL_IJ_PCG_DS | = | 'IJ-PCG-DS' | |
character(len=*), | private, | parameter | :: | HYPRE_SOL_S_PCG_NONE | = | 'PCG' |
HYPRE DOC: Preconditioned Conjugate Gradient, belongs to family of Krylov Solvers |
character(len=*), | private, | parameter | :: | HYPRE_SOL_S_PFMG_NONE | = | 'PFMG' |
HYPRE DOC: PFMG is a parallel semicoarsening multigrid solver that uses pointwise relaxation. For periodic problems, users should try to set the grid size in periodic dimensions to be as close to a power-of-two as possible. That is, if the grid size in a periodic dimension is given by N=2^m * M where M is not a power-of-two, then M should be as small as possible. Large values of will generally result in slower convergence rates. Note: PFMG is not as robust as SMG, but is much more efficient per V-cycle. |
character(len=*), | private, | parameter | :: | HYPRE_SOL_S_SMG_NONE | = | 'SMG' |
HYPRE DOC: SMG is a parallel semicoarsening multigrid solver for solving \nabla \cdot \left( D\nabla u\right) + \sigma u = f on logically rectangular grids. The code solves both 2D and 3D problems with discretization stencils of up to 9-points in 2D and up to 27-points in 3D. Note: SMG is a particularly robust method. The algorithm semicoarsens in the z-direction and uses plane smoothing. The xy plane-solves are effected by one V-cycle of the 2D SMG algorithm, which semicoarsens in the y-direction and uses line smoothing. |
A utility to call/use HYPRE scallable solvers
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
integer, | public | :: | MaxIt | = | 20 |
Maximum number of subiterations |
|
real(kind=wp), | public | :: | MaxTol | = | 1.0e-8_wp |
Maximum relative tolerance |
|
integer, | public | :: | RelaxType | = | 1 |
Relaxation Type for multigrid solvers |
|
type(block_obj), | public, | pointer | :: | block | => | null() |
Associated block structure |
integer, | public, | pointer | :: | cols(:) |
Storage array, column indices |
||
type(MPI_Comm), | public | :: | comm |
MPI communicator used by Hyper |
|||
integer, | public | :: | dim | = | 3 |
Dimensionality of the problem |
|
integer(kind=leapI8), | public | :: | grid |
Hypre grid pointer |
|||
real(kind=wp), | public | :: | h2 |
Scaling factor |
|||
integer, | public | :: | interface | = | -1 |
Hypre interface |
|
type(eulerian_obj_i), | public | :: | irow |
Row index (one row = one grid point) |
|||
integer, | public | :: | irow_hi |
low and high row indices |
|||
integer, | public | :: | irow_lo |
low and high row indices |
|||
integer, | public | :: | it | = | 0 |
Number of iteration at end of solve |
|
real(kind=wp), | public, | allocatable | :: | mat(:,:,:,:,:) |
Matrix |
||
integer(kind=leapI8), | public | :: | matrix |
Hypre matrix pointer |
|||
type(c_ptr), | public | :: | p_cols | ||||
type(c_ptr), | public | :: | p_rhs_values | ||||
type(c_ptr), | public | :: | p_rows | ||||
type(c_ptr), | public | :: | p_sol_values | ||||
type(c_ptr), | public | :: | p_tmpi | ||||
type(c_ptr), | public | :: | p_values | ||||
integer(kind=leapI8), | public | :: | par_matrix |
Pointer to ParCSR storage |
|||
integer(kind=leapI8), | public | :: | par_rhs |
Pointer to ParCSR storage |
|||
integer(kind=leapI8), | public | :: | par_sol |
Pointer to ParCSR storage |
|||
type(parallel_obj), | public, | pointer | :: | parallel | => | null() |
Associated parallel structure |
integer(kind=leapI8), | public | :: | precond |
Hypre preconditioner pointer |
|||
real(kind=wp), | public | :: | rel |
Relative error at end of solve |
|||
integer(kind=leapI8), | public | :: | rhs |
Hypre vector pointer for right-hand side |
|||
real(kind=wp), | public, | pointer | :: | rhs_values(:) |
Storage arrays, rhs |
||
integer, | public, | pointer | :: | rows(:) |
Storage arrays, row indices |
||
integer(kind=leapI8), | public | :: | sol |
Hypre vector pointer for solution |
|||
real(kind=wp), | public, | pointer | :: | sol_values(:) |
Storage arrays, sol |
||
integer(kind=leapI8), | public | :: | solver |
Hypre solver pointer |
|||
character(len=:), | public, | allocatable | :: | solver_name |
Name of solver (see list below) |
||
integer, | public | :: | st | = | 1 |
Stencil bounds (-st to +st); st=1 gives a 7-point stencil in 3D |
|
integer, | public | :: | st_size |
Stencil size |
|||
integer(kind=leapI8), | public | :: | stencil |
Hypre stencil pointer |
|||
integer, | public | :: | stm | = | 1 |
Stencil extent for matrix |
|
integer, | public, | pointer | :: | tmpi(:) |
Storage arrays, tmp integer values |
||
real(kind=wp), | public, | pointer | :: | values(:) |
Storage arrays, matrix coefficients |
procedure, public :: Finalize => hypre_obj_Final | |
procedure, public :: Initialize => hypre_obj_Init | |
procedure, public :: SelectSolver => hypre_obj_SelectSolver | |
procedure, public :: SetRHS => hypre_obj_SetRHS | |
procedure, public :: Setup => hypre_obj_Setup | |
procedure, public :: Solve => hypre_obj_Solve | |
procedure, public :: hypre_obj_BuildMatrixIJ | |
procedure, public :: hypre_obj_BuildMatrixS | |
procedure, public :: hypre_obj_PrintMatrixS | |
procedure, public :: hypre_obj_SetRHSIJ | |
procedure, public :: hypre_obj_SetRHSS | |
procedure, public :: hypre_obj_SetupGridS | |
procedure, public :: hypre_obj_SetupMatrixIJ | |
procedure, public :: hypre_obj_SetupMatrixS | |
procedure, public :: hypre_obj_SetupPointersIJ | |
procedure, public :: hypre_obj_SetupRHSIJ | |
procedure, public :: hypre_obj_SetupRHSS | |
procedure, public :: hypre_obj_SetupRowsIJ | |
procedure, public :: hypre_obj_SetupSolIJ | |
procedure, public :: hypre_obj_SetupSolIJver | |
procedure, public :: hypre_obj_SetupSolS | |
procedure, public :: hypre_obj_SetupSolverS | |
procedure, public :: hypre_obj_SetupStencilS | |
procedure, public :: hypre_obj_SolveIJ | |
procedure, public :: hypre_obj_SolveS |
Set the coefficients of the matrix
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Define the entries of the matrix Ax=b one row at a time Finite difference/Finite volume 2nd order Laplacian: ddu/dxdx = -2 u(i,j,k)/dxdx + u(i-1,j,k)/dxdx + u(i+1,j,k)/dxdx ddu/dydy = -2 u(i,j,k)/dydy + u(i,j-1,k)/dydy + u(i,j+1,k)/dydy ddu/dzdz = -2 u(i,j,k)/dzdz + u(i,j,k-1)/dzdz + u(i,j,k+1)/dzdz
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Destroy objects/pointers and clear data
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Initialize the hypre object
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(block_obj), | intent(in), | target | :: | block |
A block object |
|
type(parallel_obj), | intent(in), | target | :: | parallel |
parallel structure from main program |
Setup the hypre grid
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Select one of the preconfigured solvers and get solver-specific parameters
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
character(len=*), | intent(in) | :: | name |
Name of solver to be used |
||
real(kind=wp), | intent(in), | optional | :: | MaxTol |
Maximum relative tolerance |
|
integer, | intent(in), | optional | :: | MaxIt |
Maximum number of subiterations |
|
integer, | intent(in), | optional | :: | RelaxType |
Relaxation type |
Set the entries of the rhs vector, one element at a time
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(in) | :: | rhs |
Right hand side |
Set the entries of the rhs vector, one element at a time
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(in) | :: | rhs |
Right hand side |
Set the entries of the rhs vector, one element at a time
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(in) | :: | rhs |
Right hand side |
Setup the hypre objects in preparation for solves Note: Setting up HYPRE is an expensive operation so it's best to do this only once during a simulation as opposed to setting-up and destorying each time-step.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the hypre grid
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup matrix with IJ interface
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup and build the matrix
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the rhs vector
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the rhs vector
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup row indexing used with IJ interface
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the solution vector, and initialize it to zero
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup solver with IJ interface
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the solution vector, and initialize it to zero
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the hypre solver
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Setup the discretization stencil Each entry represents the relative offset (in index space) E.g.: a 2D 5-pt stencil would have the following geometry -- Offset { {0,0}, {-1,0}, {1,0}, {0,-1}, {0,1} } E.g.: a 3D 7-pt stencil would have the following geometry -- Offset { {0,0,0}, {-1,0,0}, {1,0,0}, {0,-1,0}, {0,1,0}, {0,0,-1}, {0,0,1} }
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
Solve the system Ax=b and return the solution
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(inout) | :: | sol |
Solution vector |
Solve the system Ax=b and return the solution, IJ interface
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(inout) | :: | sol |
Solution vector |
Solve the system Ax=b and return the solution, Struct interface
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(hypre_obj), | intent(inout) | :: | this |
Hypre machinery |
||
type(eulerian_obj_r), | intent(inout) | :: | sol |
Solution vector |