MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
PetscImsPreconditioner.F90
Go to the documentation of this file.
2 #include <petsc/finclude/petscksp.h>
3  use petscksp
4  use kindmodule, only: dp, i4b
6  use petscmatrixmodule, only: petscmatrixtype
9  use simvariablesmodule, only: iout
10 
11  implicit none
12  private
13 
14  public :: pcshellctxtype
15  public :: pcshell_apply
16  public :: pcshell_setup
17 
18  character(len=9), dimension(4) :: ims_pc_type = ["IMS ILU0 ", &
19  "IMS MILU0", &
20  "IMS ILUT ", &
21  "IMS MILUT"]
22 
24  character(len=LENMEMPATH) :: memory_path !< the memory path, taken from the system matrix
25  class(petscmatrixtype), pointer :: system_matrix !< pointer to the petsc system matrix
26  type(imslinearsettingstype), pointer :: linear_settings !< the linear settings from IMS
27  integer(I4B) :: ipc !< the IMS preconditioner type: 1=ILU0, 2=MILU0, 3=ILUT, 4=MILUT
28  character(len=16) :: ims_pc_type !< the IMS preconditioner type as a string
29  ! PC matrix
30  integer(I4B), dimension(:), contiguous, pointer :: iapc => null() !< CSR row pointer
31  integer(I4B), dimension(:), contiguous, pointer :: japc => null() !< CSR columns
32  real(dp), dimension(:), contiguous, pointer :: apc => null() !< CSR Preconditioner matrix values
33  ! work vectors
34  integer(I4B), dimension(:), contiguous, pointer :: iw => null() !< (M)ILU0 work array
35  real(dp), dimension(:), contiguous, pointer :: w => null() !< (M)ILU0 work array
36  integer(I4B), dimension(:), contiguous, pointer :: jlu => null() !< (M)ILUT work array
37  integer(I4B), dimension(:), contiguous, pointer :: jw => null() !< (M)ILUT work array
38  real(dp), dimension(:), contiguous, pointer :: wlu => null() !< (M)ILUT work array
39  contains
40  procedure :: create => pctx_create
41  procedure :: destroy => pctx_destroy
42 
43  end type
44 
45  interface
46  subroutine pcshellgetcontext(pc, ctx, ierr)
47  import pcshellctxtype, tpc
48  type(tpc) :: pc
49  type(pcshellctxtype), pointer :: ctx
50  integer :: ierr
51  end subroutine
52  end interface
53 
54 contains
55 
56  !> @brief Create the preconditioner context from the system matrix
57  !<
58  subroutine pctx_create(this, matrix, lin_settings)
61  class(pcshellctxtype) :: this !< this instance
62  class(petscmatrixtype), pointer :: matrix !< the linear system matrix
63  type(imslinearsettingstype), pointer :: lin_settings !< linear settings from IMS
64  ! local
65  integer(I4B) :: n, neq, nja
66  integer(I4B), dimension(:), contiguous, pointer :: ia, ja
67  real(DP), dimension(:), contiguous, pointer :: amat
68  integer(I4B) :: niapc, njapc, njlu, njw, nwlu
69 
70  this%memory_path = matrix%memory_path
71  this%linear_settings => lin_settings
72 
73  this%ipc = 1
74  if (this%linear_settings%level > 0) this%ipc = this%ipc + 2
75  if (this%linear_settings%relax > 0.0_dp) this%ipc = this%ipc + 1
76  this%ims_pc_type = ims_pc_type(this%ipc)
77 
78  this%system_matrix => matrix
79  call matrix%get_aij_local(ia, ja, amat)
80 
81  neq = matrix%nrow
82  nja = matrix%nnz_local
83  call ims_calc_pcdims(neq, nja, ia, this%linear_settings%level, this%ipc, &
84  niapc, njapc, njlu, njw, nwlu)
85 
86  ! preconditioner matrix
87  call mem_allocate(this%IAPC, niapc + 1, 'IAPC', this%memory_path)
88  call mem_allocate(this%JAPC, njapc, 'JAPC', this%memory_path)
89  call mem_allocate(this%APC, njapc, 'APC', this%memory_path)
90  if (this%ipc == 1 .or. this%ipc == 2) then
91  call ims_base_pccrs(niapc, njapc, ia, ja, this%IAPC, this%JAPC)
92  end if
93  do n = 1, njapc
94  this%APC(n) = dzero
95  end do
96 
97  ! (M)ILU0 work arrays
98  call mem_allocate(this%IW, niapc, 'IW', this%memory_path)
99  call mem_allocate(this%W, niapc, 'W', this%memory_path)
100  do n = 1, niapc
101  this%IW(n) = izero
102  this%W(n) = dzero
103  end do
104 
105  ! (M)ILUT work arrays
106  call mem_allocate(this%JLU, njlu, 'JLU', this%memory_path)
107  call mem_allocate(this%JW, njw, 'JW', this%memory_path)
108  call mem_allocate(this%WLU, nwlu, 'WLU', this%memory_path)
109 
110  do n = 1, njlu
111  this%JLU(n) = izero
112  end do
113  do n = 1, njw
114  this%JW(n) = izero
115  end do
116  do n = 1, nwlu
117  this%WLU(n) = dzero
118  end do
119 
120  end subroutine pctx_create
121 
122  !> @brief Cleanup the context
123  !<
124  subroutine pctx_destroy(this)
126  class(pcshellctxtype) :: this
127 
128  ! matrix
129  call mem_deallocate(this%IAPC)
130  call mem_deallocate(this%JAPC)
131  call mem_deallocate(this%APC)
132 
133  ! work arrays
134  call mem_deallocate(this%IW)
135  call mem_deallocate(this%W)
136  call mem_deallocate(this%JLU)
137  call mem_deallocate(this%JW)
138  call mem_deallocate(this%WLU)
139 
140  this%system_matrix => null()
141 
142  end subroutine pctx_destroy
143 
144  !> @brief Apply shell preconditioner
145  !< (this is not a type bound procedure)
146  subroutine pcshell_apply(pc, x, y, ierr)
148  external lusol ! from ilut.f90
149  pc :: pc !< the shell preconditioner
150  vec :: x !< the input vector
151  vec :: y !< the output vector
152  petscerrorcode :: ierr !< PETSc error code
153  ! local
154  type(pcshellctxtype), pointer :: pc_ctx => null()
155  real(dp), dimension(:), pointer :: local_x, local_y
156  integer(I4B) :: neq, nja
157 
158  call pcshellgetcontext(pc, pc_ctx, ierr)
159  chkerrq(ierr)
160 
161  call vecgetarrayreadf90(x, local_x, ierr)
162  chkerrq(ierr)
163  call vecgetarrayf90(y, local_y, ierr)
164  chkerrq(ierr)
165 
166  neq = pc_ctx%system_matrix%nrow
167  nja = pc_ctx%system_matrix%nnz_local
168  select case (pc_ctx%ipc)
169  case (1, 2)
170  call ims_base_ilu0a(nja, neq, pc_ctx%APC, pc_ctx%IAPC, pc_ctx%JAPC, &
171  local_x, local_y)
172  case (3, 4)
173  call lusol(neq, local_x, local_y, pc_ctx%APC, pc_ctx%JLU, pc_ctx%IW)
174  end select
175 
176  call vecrestorearrayf90(x, local_x, ierr)
177  chkerrq(ierr)
178  call vecrestorearrayf90(y, local_y, ierr)
179  chkerrq(ierr)
180 
181  end subroutine pcshell_apply
182 
183  !> @brief Set up the custom preconditioner
184  !< (this is not a type bound procedure)
185  subroutine pcshell_setup(pc, ierr)
187  pc :: pc !< the shell preconditioner
188  petscerrorcode :: ierr !< PETSc error code
189  ! local
190  type(pcshellctxtype), pointer :: pc_ctx => null()
191  integer(I4B) :: neq, nja
192  integer(I4B) :: niapc, njapc, njlu, njw, nwlu
193  integer(I4B), dimension(:), contiguous, pointer :: ia, ja
194  real(dp), dimension(:), contiguous, pointer :: amat
195 
196  call pcshellgetcontext(pc, pc_ctx, ierr)
197  chkerrq(ierr)
198 
199  ! note the two different matrix types here:
200  ! bridging between PETSc and IMS
201  call pc_ctx%system_matrix%get_aij_local(ia, ja, amat)
202 
203  neq = pc_ctx%system_matrix%nrow
204  nja = pc_ctx%system_matrix%nnz_local
205  niapc = size(pc_ctx%IAPC) - 1
206  njapc = size(pc_ctx%JAPC)
207  njlu = size(pc_ctx%JLU)
208  njw = size(pc_ctx%JW)
209  nwlu = size(pc_ctx%WLU)
210  call ims_base_pcu(iout, nja, neq, niapc, njapc, &
211  pc_ctx%ipc, pc_ctx%linear_settings%relax, &
212  amat, ia, ja, &
213  pc_ctx%APC, pc_ctx%IAPC, pc_ctx%JAPC, &
214  pc_ctx%IW, pc_ctx%W, &
215  pc_ctx%linear_settings%level, &
216  pc_ctx%linear_settings%droptol, &
217  njlu, njw, nwlu, &
218  pc_ctx%JLU, pc_ctx%JW, pc_ctx%WLU)
219 
220  end subroutine pcshell_setup
221 
subroutine lusol(n, y, x, alu, jlu, ju)
Definition: ilut.f90:466
This module contains simulation constants.
Definition: Constants.f90:9
integer(i4b), parameter izero
integer constant zero
Definition: Constants.f90:50
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64
integer(i4b), parameter lenmempath
maximum length of the memory path
Definition: Constants.f90:26
This module contains the IMS linear accelerator subroutines.
subroutine ims_base_pcu(IOUT, NJA, NEQ, NIAPC, NJAPC, IPC, RELAX, AMAT, IA, JA, APC, IAPC, JAPC, IW, W, LEVEL, DROPTOL, NJLU, NJW, NWLU, JLU, JW, WLU)
@ brief Update the preconditioner
subroutine ims_base_pccrs(NEQ, NJA, IA, JA, IAPC, JAPC)
@ brief Generate CRS pointers for the preconditioner
subroutine ims_calc_pcdims(neq, nja, ia, level, ipc, niapc, njapc, njlu, njw, nwlu)
subroutine ims_base_ilu0a(NJA, NEQ, APC, IAPC, JAPC, R, D)
@ brief Apply the ILU0 and MILU0 preconditioners
This module defines variable data types.
Definition: kind.f90:8
subroutine, public pcshell_setup(pc, ierr)
Set up the custom preconditioner.
subroutine pctx_destroy(this)
Cleanup the context.
subroutine, public pcshell_apply(pc, x, y, ierr)
Apply shell preconditioner.
subroutine pctx_create(this, matrix, lin_settings)
Create the preconditioner context from the system matrix.
character(len=9), dimension(4) ims_pc_type
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iout
file unit number for simulation output