MODFLOW 6  version 6.6.0.dev0
USGS Modular Hydrologic Model
petscimspreconditionermodule Module Reference

Data Types

type  pcshellctxtype
 
interface  PCShellGetContext
 

Functions/Subroutines

subroutine pctx_create (this, matrix, lin_settings)
 Create the preconditioner context from the system matrix. More...
 
subroutine pctx_destroy (this)
 Cleanup the context. More...
 
subroutine, public pcshell_apply (pc, x, y, ierr)
 Apply shell preconditioner. More...
 
subroutine, public pcshell_setup (pc, ierr)
 Set up the custom preconditioner. More...
 

Variables

character(len=9), dimension(4) ims_pc_type = ["IMS ILU0 ", "IMS MILU0", "IMS ILUT ", "IMS MILUT"]
 

Function/Subroutine Documentation

◆ pcshell_apply()

subroutine, public petscimspreconditionermodule::pcshell_apply (   pc,
  x,
  y,
  ierr 
)

Definition at line 146 of file PetscImsPreconditioner.F90.

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 
subroutine lusol(n, y, x, alu, jlu, ju)
Definition: ilut.f90:466
This module contains the IMS linear accelerator subroutines.
subroutine ims_base_ilu0a(NJA, NEQ, APC, IAPC, JAPC, R, D)
@ brief Apply the ILU0 and MILU0 preconditioners
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pcshell_setup()

subroutine, public petscimspreconditionermodule::pcshell_setup (   pc,
  ierr 
)

Definition at line 185 of file PetscImsPreconditioner.F90.

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 
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
Here is the call graph for this function:
Here is the caller graph for this function:

◆ pctx_create()

subroutine petscimspreconditionermodule::pctx_create ( class(pcshellctxtype this,
class(petscmatrixtype), pointer  matrix,
type(imslinearsettingstype), pointer  lin_settings 
)
private
Parameters
thisthis instance
matrixthe linear system matrix
lin_settingslinear settings from IMS

Definition at line 58 of file PetscImsPreconditioner.F90.

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 
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)
Here is the call graph for this function:

◆ pctx_destroy()

subroutine petscimspreconditionermodule::pctx_destroy ( class(pcshellctxtype this)

Definition at line 124 of file PetscImsPreconditioner.F90.

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 

Variable Documentation

◆ ims_pc_type

character(len=9), dimension(4) petscimspreconditionermodule::ims_pc_type = ["IMS ILU0 ", "IMS MILU0", "IMS ILUT ", "IMS MILUT"]
private

Definition at line 18 of file PetscImsPreconditioner.F90.

18  character(len=9), dimension(4) :: IMS_PC_TYPE = ["IMS ILU0 ", &
19  "IMS MILU0", &
20  "IMS ILUT ", &
21  "IMS MILUT"]