hypre_obj Derived Type

type, public :: hypre_obj

A utility to call/use HYPRE scallable solvers


Inherits

type~~hypre_obj~~InheritsGraph type~hypre_obj hypre_obj MPI_Comm MPI_Comm type~hypre_obj->MPI_Comm comm c_ptr c_ptr type~hypre_obj->c_ptr p_rhs_values, p_sol_values, p_values, p_rows, p_cols, p_tmpi type~block_obj block_obj type~hypre_obj->type~block_obj block type~eulerian_obj_i eulerian_obj_i type~hypre_obj->type~eulerian_obj_i irow type~parallel_obj parallel_obj type~hypre_obj->type~parallel_obj parallel type~block_obj->type~parallel_obj parallel MPI_Datatype MPI_Datatype type~block_obj->MPI_Datatype gc_slab_r, gc_slab_i type~axis_obj axis_obj type~block_obj->type~axis_obj axis, axis_partition type~eulerian_obj_base eulerian_obj_base type~eulerian_obj_i->type~eulerian_obj_base type~parallel_obj->MPI_Datatype REAL_SP, REAL_DP, REAL_WP, COMPLEX_SP, COMPLEX_DP, COMPLEX_WP, INTEGER, INT8, LOGICAL MPI_Info MPI_Info type~parallel_obj->MPI_Info mpi_info type~communicators communicators type~parallel_obj->type~communicators comm type~patch patch type~parallel_obj->type~patch rank type~communicators->MPI_Comm w, g type~eulerian_obj_base->type~block_obj block type~eulerian_obj_base->type~parallel_obj parallel

Inherited by

type~~hypre_obj~~InheritedByGraph type~hypre_obj hypre_obj type~cdifs_obj cdifs_obj type~cdifs_obj->type~hypre_obj hypre type~grans_obj grans_obj type~grans_obj->type~hypre_obj hypre

Components

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


Type-Bound Procedures

procedure, public :: Finalize => hypre_obj_Final

  • private subroutine hypre_obj_Final(this)

    Destroy objects/pointers and clear data

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: Initialize => hypre_obj_Init

  • private subroutine hypre_obj_Init(this, block, parallel)

    Initialize the hypre object

    Arguments

    Type IntentOptional 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

procedure, public :: SelectSolver => hypre_obj_SelectSolver

  • private subroutine hypre_obj_SelectSolver(this, name, MaxTol, MaxIt, RelaxType)

    Select one of the preconfigured solvers and get solver-specific parameters

    Arguments

    Type IntentOptional 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

procedure, public :: SetRHS => hypre_obj_SetRHS

  • private subroutine hypre_obj_SetRHS(this, rhs)

    Set the entries of the rhs vector, one element at a time

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(in) :: rhs

    Right hand side

procedure, public :: Setup => hypre_obj_Setup

  • private subroutine hypre_obj_Setup(this)

    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.

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: Solve => hypre_obj_Solve

  • private subroutine hypre_obj_Solve(this, sol)

    Solve the system Ax=b and return the solution

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(inout) :: sol

    Solution vector

procedure, public :: hypre_obj_BuildMatrixIJ

  • private subroutine hypre_obj_BuildMatrixIJ(this)

    Set the coefficients of the matrix

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_BuildMatrixS

  • private subroutine hypre_obj_BuildMatrixS(this)

    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

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_PrintMatrixS

  • private subroutine hypre_obj_PrintMatrixS(this)

    Setup the hypre grid

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetRHSIJ

  • private subroutine hypre_obj_SetRHSIJ(this, rhs)

    Set the entries of the rhs vector, one element at a time

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(in) :: rhs

    Right hand side

procedure, public :: hypre_obj_SetRHSS

  • private subroutine hypre_obj_SetRHSS(this, rhs)

    Set the entries of the rhs vector, one element at a time

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(in) :: rhs

    Right hand side

procedure, public :: hypre_obj_SetupGridS

  • private subroutine hypre_obj_SetupGridS(this)

    Setup the hypre grid

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupMatrixIJ

  • private subroutine hypre_obj_SetupMatrixIJ(this)

    Setup matrix with IJ interface

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupMatrixS

procedure, public :: hypre_obj_SetupPointersIJ

procedure, public :: hypre_obj_SetupRHSIJ

  • private subroutine hypre_obj_SetupRHSIJ(this)

    Setup the rhs vector

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupRHSS

  • private subroutine hypre_obj_SetupRHSS(this)

    Setup the rhs vector

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupRowsIJ

  • private subroutine hypre_obj_SetupRowsIJ(this)

    Setup row indexing used with IJ interface

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupSolIJ

  • private subroutine hypre_obj_SetupSolIJ(this)

    Setup the solution vector, and initialize it to zero

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupSolIJver

  • private subroutine hypre_obj_SetupSolIJver(this)

    Setup solver with IJ interface

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupSolS

  • private subroutine hypre_obj_SetupSolS(this)

    Setup the solution vector, and initialize it to zero

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupSolverS

  • private subroutine hypre_obj_SetupSolverS(this)

    Setup the hypre solver

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SetupStencilS

  • private subroutine hypre_obj_SetupStencilS(this)

    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} }

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

procedure, public :: hypre_obj_SolveIJ

  • private subroutine hypre_obj_SolveIJ(this, sol)

    Solve the system Ax=b and return the solution, IJ interface

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(inout) :: sol

    Solution vector

procedure, public :: hypre_obj_SolveS

  • private subroutine hypre_obj_SolveS(this, sol)

    Solve the system Ax=b and return the solution, Struct interface

    Arguments

    Type IntentOptional Attributes Name
    class(hypre_obj), intent(inout) :: this

    Hypre machinery

    type(eulerian_obj_r), intent(inout) :: sol

    Solution vector