MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
LinearSolverBase.f90
Go to the documentation of this file.
2  use kindmodule, only: i4b, dp
8 
9  implicit none
10  private
11 
12  !> @brief Abstract type for linear solver
13  !!
14  !! This serves as the base type for our solvers:
15  !! sequential, parallel, petsc, block solver, ...
16  !<
17  type, public, abstract :: linearsolverbasetype
18  character(len=LENSOLUTIONNAME) :: name
19  integer(I4B) :: nitermax
20  integer(I4B) :: iteration_number
21  integer(I4B) :: is_converged
22  contains
23  procedure(initialize_if), deferred :: initialize
24  procedure(print_summary_if), deferred :: print_summary
25  procedure(solve_if), deferred :: solve
26  procedure(destroy_if), deferred :: destroy
27 
28  procedure(create_matrix_if), deferred :: create_matrix
29  end type linearsolverbasetype
30 
31  abstract interface
32  subroutine initialize_if(this, matrix, linear_settings, convergence_summary)
35  class(linearsolverbasetype) :: this
36  class(matrixbasetype), pointer :: matrix
37  type(imslinearsettingstype), pointer :: linear_settings
38  type(convergencesummarytype), pointer :: convergence_summary
39  end subroutine
40  subroutine print_summary_if(this)
42  class(linearsolverbasetype) :: this
43  end subroutine
44  subroutine solve_if(this, kiter, rhs, x, cnvg_summary)
46  class(linearsolverbasetype) :: this
47  integer(I4B) :: kiter
48  class(vectorbasetype), pointer :: rhs
49  class(vectorbasetype), pointer :: x
50  type(convergencesummarytype) :: cnvg_summary
51  end subroutine
52  subroutine destroy_if(this)
54  class(linearsolverbasetype) :: this
55  end subroutine
56  function create_matrix_if(this) result(matrix)
58  class(linearsolverbasetype) :: this
59  class(matrixbasetype), pointer :: matrix
60  end function
61  end interface
62 
63 end module linearsolverbasemodule
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter lensolutionname
maximum length of the solution name
Definition: Constants.f90:21
subroutine destroy(this)
Cleanup.
This module defines variable data types.
Definition: kind.f90:8
This structure stores the generic convergence info for a solution.