petsc_controller_class Module

Houses class in charge of interfacing with PETSc


Used by


Interfaces

interface

  • private pure module function objectsCreated(this) result(created)

    Return true if PETSc objects have been created

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(in) :: this

    Return Value logical

interface

  • private pure module function getLastConvergedReason(this) result(conv)

    Getter for lastConvergedReason

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(in) :: this

    Return Value integer(kind=ik)

interface

  • private pure module function getLastNumIterations(this) result(numIters)

    Getter for lastNumIterations

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(in) :: this

    Return Value integer(kind=ik)

interface

  • private module subroutine initPETScController(this, indexingObj, mpiCont, solOptions, numObjs)

    PETSc controller initialization routine

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(inout) :: this
    type(Indexing), intent(in) :: indexingObj

    Indexing object used to initialize preallocation data

    type(MPIController), intent(in) :: mpiCont

    Reference MPI controller

    type(SolverOptions), intent(in) :: solOptions

    KSP solver options

    integer(kind=ik), intent(in), optional :: numObjs

    Number of PETSc object groups

interface

  • private pure module subroutine addRowDataToPreallocation(this, rowData)

    Add sparse row data structure to total preallocation sparsity pattern

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(inout) :: this
    type(SparseRowData), intent(in) :: rowData

interface

  • private module subroutine addRowValuesToMatrix(this, rowData, multConst, objGroup)

    Add row data to matrix, optionally multiplied by a real constant

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(inout) :: this
    type(SparseRowData), intent(in) :: rowData
    real(kind=rk), intent(in), optional :: multConst
    integer(kind=ik), intent(in), optional :: objGroup

interface

  • private module subroutine createPETScObjs(this)

    Assemble the preallocation object, create PETSc objects, and deallocate the preallocation object

    Arguments

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

interface

  • private module subroutine linearSolve(this, knownVec, unknownVec, objGroup)

    Solve the equation petscMat * unknownVec = knownVec where the vectors are local to each processor - assembles and zeros out the matrix before/after the solve

    Arguments

    Type IntentOptional Attributes Name
    class(PETScController), intent(inout) :: this
    real(kind=rk), intent(in), dimension(:) :: knownVec
    real(kind=rk), intent(out), dimension(:) :: unknownVec
    integer(kind=ik), intent(in), optional :: objGroup

interface

  • private module subroutine finalize(this)

    Finalizes PETSc and destroys all PETSc objects

    Arguments

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

Derived Types

type, public ::  SolverOptions

Contains data used by the PETSc ksp solver to set various options (refer to PETSc manual)

Components

Type Visibility Attributes Name Initial
character(len=:), public, allocatable :: kspSolverType

String representing the ksp solver type

real(kind=rk), public :: solverToleranceRel = real(1.0d-17, kind=rk)

Relative solver tolerance

real(kind=rk), public :: solverToleranceAbs = real(1.0d-20, kind=rk)

Absolute solver tolerance

real(kind=rk), public :: solverToleranceDiv = real(1.0d07, kind=rk)

Solver divergence tolerance

integer(kind=ik), public :: maxSolverIters = 10000

Maximum number of solver iterations

character(len=:), public, allocatable :: hyprePC

String representing the hypre PC type

character(len=:), public, allocatable :: petscOptions

String representing command line style options

type, public ::  PETSsObjs

Components

Type Visibility Attributes Name Initial
vec, private :: rhs

RHS PETSc vector object in linear solve

vec, private :: sol

LHS PETSc vector object in linear solve

mat, private :: petscMat

PETSc matrix object used in linear solve

pc, private :: preconditioner

Precondition object used in solve

ksp, private :: solver

KSP solver object

type, public, extends(Object) ::  PETScController

Object responsible for interfacing with PETSc

Components

Type Visibility Attributes Name Initial
logical, public :: userDefined = .false.

True only if user explicitly sets it to true

type(PETScPreallocationData), private :: preallocData

Data used for preallocating PETSc objects

type(SolverOptions), private :: options

Solver options for the KSP solver

logical, private :: objCreated

True if PETSc object components have been created

type(PETSsObjs), private, allocatable, dimension(:) :: objs

Groups of PETSc objects

ksp, public :: ConvergedReason::lastConvergedReason

Converged reason for last attempted solve

petscint, public :: lastNumIterations

Number of iterations for last attemtped solve

Type-Bound Procedures

procedure, public :: isDefined => isDefinedObject
procedure, public :: makeDefined => makeDefinedObject
procedure, public :: makeUndefined => makeUndefinedObject
procedure, public :: addRowDataToPreallocation
procedure, public :: createPETScObjs
procedure, public :: addRowValuesToMatrix
procedure, public :: linearSolve
procedure, public :: objectsCreated
procedure, public :: getLastConvergedReason
procedure, public :: getLastNumIterations
procedure, public :: init => initPETScController
procedure, public :: finalize