MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
imslinearmodule Module Reference

Data Types

type  imslineardatatype
 

Functions/Subroutines

subroutine imslinear_ar (this, NAME, IOUT, IPRIMS, MXITER, NEQ, matrix, RHS, X, linear_settings)
 @ brief Allocate storage and read data More...
 
subroutine imslinear_summary (this, mxiter)
 @ brief Write summary of settings More...
 
subroutine allocate_scalars (this)
 @ brief Allocate and initialize scalars More...
 
subroutine imslinear_da (this)
 @ brief Deallocate memory More...
 
subroutine imslinear_set_input (this, IFDPARAM)
 @ brief Set default settings More...
 
subroutine imslinear_ap (this, ICNVG, KSTP, KITER, IN_ITER, NCONV, CONVNMOD, CONVMODSTART, CACCEL, summary)
 @ brief Base linear accelerator subroutine More...
 

Function/Subroutine Documentation

◆ allocate_scalars()

subroutine imslinearmodule::allocate_scalars ( class(imslineardatatype), intent(inout)  this)
private

Allocate and initialize linear accelerator scalars

Parameters
[in,out]thisImsLinearDataType instance

Definition at line 455 of file ImsLinear.f90.

456  ! -- modules
458  ! -- dummy variables
459  class(ImsLinearDataType), intent(inout) :: this !< ImsLinearDataType instance
460  !
461  ! -- allocate scalars
462  call mem_allocate(this%iout, 'IOUT', this%memoryPath)
463  call mem_allocate(this%ipc, 'IPC', this%memoryPath)
464  call mem_allocate(this%iacpc, 'IACPC', this%memoryPath)
465  call mem_allocate(this%niterc, 'NITERC', this%memoryPath)
466  call mem_allocate(this%niabcgs, 'NIABCGS', this%memoryPath)
467  call mem_allocate(this%niapc, 'NIAPC', this%memoryPath)
468  call mem_allocate(this%njapc, 'NJAPC', this%memoryPath)
469  call mem_allocate(this%epfact, 'EPFACT', this%memoryPath)
470  call mem_allocate(this%l2norm0, 'L2NORM0', this%memoryPath)
471  call mem_allocate(this%njlu, 'NJLU', this%memoryPath)
472  call mem_allocate(this%njw, 'NJW', this%memoryPath)
473  call mem_allocate(this%nwlu, 'NWLU', this%memoryPath)
474  !
475  ! -- initialize scalars
476  this%iout = 0
477  this%ipc = 0
478  this%iacpc = 0
479  this%niterc = 0
480  this%niabcgs = 0
481  this%niapc = 0
482  this%njapc = 0
483  this%epfact = dzero
484  this%l2norm0 = 0
485  this%njlu = 0
486  this%njw = 0
487  this%nwlu = 0

◆ imslinear_ap()

subroutine imslinearmodule::imslinear_ap ( class(imslineardatatype), intent(inout)  this,
integer(i4b), intent(inout)  ICNVG,
integer(i4b), intent(in)  KSTP,
integer(i4b), intent(in)  KITER,
integer(i4b), intent(inout)  IN_ITER,
integer(i4b), intent(in)  NCONV,
integer(i4b), intent(in)  CONVNMOD,
integer(i4b), dimension(convnmod + 1), intent(inout)  CONVMODSTART,
character(len=31), dimension(nconv), intent(inout)  CACCEL,
type(convergencesummarytype), intent(in), pointer  summary 
)
private

Base linear accelerator subroutine that scales and reorders the system of equations, if necessary, updates the preconditioner, and calls the appropriate linear accelerator.

Parameters
[in,out]thisImsLinearDataType instance
[in,out]icnvgconvergence flag (1) non-convergence (0)
[in]kstptime step number
[in]kiterouter iteration number
[in,out]in_iterinner iteration number
[in]summaryConvergence summary report

Definition at line 617 of file ImsLinear.f90.

620  ! -- modules
621  USE simmodule
622  ! -- dummy variables
623  CLASS(ImsLinearDataType), INTENT(INOUT) :: this !< ImsLinearDataType instance
624  integer(I4B), INTENT(INOUT) :: ICNVG !< convergence flag (1) non-convergence (0)
625  integer(I4B), INTENT(IN) :: KSTP !< time step number
626  integer(I4B), INTENT(IN) :: KITER !< outer iteration number
627  integer(I4B), INTENT(INOUT) :: IN_ITER !< inner iteration number
628  ! -- convergence information dummy variables
629  integer(I4B), INTENT(IN) :: NCONV !<
630  integer(I4B), INTENT(IN) :: CONVNMOD !<
631  integer(I4B), DIMENSION(CONVNMOD + 1), INTENT(INOUT) :: CONVMODSTART !<
632  character(len=31), DIMENSION(NCONV), INTENT(INOUT) :: CACCEL !<
633  type(ConvergenceSummaryType), pointer, intent(in) :: summary !< Convergence summary report
634  ! -- local variables
635  integer(I4B) :: n
636  integer(I4B) :: innerit
637  integer(I4B) :: irc
638  integer(I4B) :: itmax
639  real(DP) :: dnrm2
640  !
641  ! -- set epfact based on timestep
642  this%EPFACT = ims_base_epfact(this%ICNVGOPT, kstp)
643  !
644  ! -- SCALE PROBLEM
645  IF (this%ISCL .NE. 0) THEN
646  CALL ims_base_scale(0, this%ISCL, &
647  this%NEQ, this%NJA, this%IA, this%JA, &
648  this%AMAT, this%X, this%RHS, &
649  this%DSCALE, this%DSCALE2)
650  END IF
651  !
652  ! -- PERMUTE ROWS, COLUMNS, AND RHS
653  IF (this%IORD /= 0) THEN
654  CALL dperm(this%NEQ, this%AMAT, this%JA, this%IA, &
655  this%ARO, this%JARO, this%IARO, &
656  this%LORDER, this%ID, 1)
657  CALL dvperm(this%NEQ, this%X, this%LORDER)
658  CALL dvperm(this%NEQ, this%RHS, this%LORDER)
659  this%IA0 => this%IARO
660  this%JA0 => this%JARO
661  this%A0 => this%ARO
662  ELSE
663  this%IA0 => this%IA
664  this%JA0 => this%JA
665  this%A0 => this%AMAT
666  END IF
667  !
668  ! -- UPDATE PRECONDITIONER
669  CALL ims_base_pcu(this%iout, this%NJA, this%NEQ, this%NIAPC, this%NJAPC, &
670  this%IPC, this%RELAX, this%A0, this%IA0, this%JA0, &
671  this%APC, this%IAPC, this%JAPC, this%IW, this%W, &
672  this%LEVEL, this%DROPTOL, this%NJLU, this%NJW, &
673  this%NWLU, this%JLU, this%JW, this%WLU)
674  !
675  ! -- INITIALIZE SOLUTION VARIABLE AND ARRAYS
676  IF (kiter == 1) then
677  this%NITERC = 0
678  summary%iter_cnt = 0
679  end if
680  irc = 1
681  icnvg = 0
682  DO n = 1, this%NEQ
683  this%D(n) = dzero
684  this%P(n) = dzero
685  this%Q(n) = dzero
686  this%Z(n) = dzero
687  END DO
688  !
689  ! -- CALCULATE INITIAL RESIDUAL
690  call ims_base_residual(this%NEQ, this%NJA, this%X, this%RHS, this%D, &
691  this%A0, this%IA0, this%JA0)
692  this%L2NORM0 = dnrm2(this%NEQ, this%D, 1)
693  !
694  ! -- CHECK FOR EXACT SOLUTION
695  itmax = this%ITER1
696  IF (this%L2NORM0 == dzero) THEN
697  itmax = 0
698  icnvg = 1
699  END IF
700  !
701  ! -- SOLUTION BY THE CONJUGATE GRADIENT METHOD
702  IF (this%ILINMETH == 1) THEN
703  CALL ims_base_cg(icnvg, itmax, innerit, &
704  this%NEQ, this%NJA, this%NIAPC, this%NJAPC, &
705  this%IPC, this%ICNVGOPT, this%NORTH, &
706  this%DVCLOSE, this%RCLOSE, this%L2NORM0, &
707  this%EPFACT, this%IA0, this%JA0, this%A0, &
708  this%IAPC, this%JAPC, this%APC, &
709  this%X, this%RHS, this%D, this%P, this%Q, this%Z, &
710  this%NJLU, this%IW, this%JLU, &
711  nconv, convnmod, convmodstart, &
712  caccel, summary)
713  !
714  ! -- SOLUTION BY THE BICONJUGATE GRADIENT STABILIZED METHOD
715  ELSE IF (this%ILINMETH == 2) THEN
716  CALL ims_base_bcgs(icnvg, itmax, innerit, &
717  this%NEQ, this%NJA, this%NIAPC, this%NJAPC, &
718  this%IPC, this%ICNVGOPT, this%NORTH, &
719  this%ISCL, this%DSCALE, &
720  this%DVCLOSE, this%RCLOSE, this%L2NORM0, &
721  this%EPFACT, this%IA0, this%JA0, this%A0, &
722  this%IAPC, this%JAPC, this%APC, &
723  this%X, this%RHS, this%D, this%P, this%Q, &
724  this%T, this%V, this%DHAT, this%PHAT, this%QHAT, &
725  this%NJLU, this%IW, this%JLU, &
726  nconv, convnmod, convmodstart, &
727  caccel, summary)
728  END IF
729  !
730  ! -- BACK PERMUTE AMAT, SOLUTION, AND RHS
731  IF (this%IORD /= 0) THEN
732  CALL dperm(this%NEQ, this%A0, this%JA0, this%IA0, &
733  this%AMAT, this%JA, this%IA, &
734  this%IORDER, this%ID, 1)
735  CALL dvperm(this%NEQ, this%X, this%IORDER)
736  CALL dvperm(this%NEQ, this%RHS, this%IORDER)
737  END IF
738  !
739  ! -- UNSCALE PROBLEM
740  IF (this%ISCL .NE. 0) THEN
741  CALL ims_base_scale(1, this%ISCL, &
742  this%NEQ, this%NJA, this%IA, this%JA, &
743  this%AMAT, this%X, this%RHS, &
744  this%DSCALE, this%DSCALE2)
745  END IF
746  !
747  ! -- SET IMS INNER ITERATION NUMBER (IN_ITER) TO NUMBER OF
748  ! IMSLINEAR INNER ITERATIONS (innerit)
749  in_iter = innerit
real(kind=8) function dnrm2(n, x, incx)
Definition: blas1_d.f90:388
This module contains simulation methods.
Definition: Sim.f90:10
subroutine dvperm(n, x, perm)
Definition: sparsekit.f90:62
subroutine dperm(nrow, a, ja, ia, ao, jao, iao, perm, qperm, job)
Definition: sparsekit.f90:354
Here is the call graph for this function:

◆ imslinear_ar()

subroutine imslinearmodule::imslinear_ar ( class(imslineardatatype), intent(inout)  this,
character(len=lensolutionname), intent(in)  NAME,
integer(i4b), intent(in)  IOUT,
integer(i4b), intent(in), target  IPRIMS,
integer(i4b), intent(in)  MXITER,
integer(i4b), intent(in), target  NEQ,
class(matrixbasetype), pointer  matrix,
real(dp), dimension(neq), intent(inout), target  RHS,
real(dp), dimension(neq), intent(inout), target  X,
type(imslinearsettingstype), pointer  linear_settings 
)

Allocate storage for linear accelerators and read data

Parameters
[in,out]thisImsLinearDataType instance
[in]namesolution name
[in]ioutsimulation listing file unit
[in]iprimsprint option
[in]mxitermaximum outer iterations
[in]neqnumber of equations
[in,out]rhsright-hand side
[in,out]xdependent variables
linear_settingsthe settings form the IMS file

Definition at line 111 of file ImsLinear.f90.

113  ! -- modules
116  use simmodule, only: store_error, count_errors, &
118  ! -- dummy variables
119  CLASS(ImsLinearDataType), INTENT(INOUT) :: this !< ImsLinearDataType instance
120  CHARACTER(LEN=LENSOLUTIONNAME), INTENT(IN) :: NAME !< solution name
121  integer(I4B), INTENT(IN) :: IOUT !< simulation listing file unit
122  integer(I4B), TARGET, INTENT(IN) :: IPRIMS !< print option
123  integer(I4B), INTENT(IN) :: MXITER !< maximum outer iterations
124  integer(I4B), TARGET, INTENT(IN) :: NEQ !< number of equations
125  class(MatrixBaseType), pointer :: matrix
126  real(DP), DIMENSION(NEQ), TARGET, INTENT(INOUT) :: RHS !< right-hand side
127  real(DP), DIMENSION(NEQ), TARGET, INTENT(INOUT) :: X !< dependent variables
128  type(ImsLinearSettingsType), pointer :: linear_settings !< the settings form the IMS file
129  ! -- local variables
130  character(len=LINELENGTH) :: errmsg
131  integer(I4B) :: n
132  integer(I4B) :: i0
133  integer(I4B) :: iscllen, iolen
134 
135  !
136  ! -- DEFINE NAME
137  this%memoryPath = create_mem_path(name, 'IMSLINEAR')
138  !
139  ! -- SET pointers to IMS settings
140  this%DVCLOSE => linear_settings%dvclose
141  this%RCLOSE => linear_settings%rclose
142  this%ICNVGOPT => linear_settings%icnvgopt
143  this%ITER1 => linear_settings%iter1
144  this%ILINMETH => linear_settings%ilinmeth
145  this%iSCL => linear_settings%iscl
146  this%IORD => linear_settings%iord
147  this%NORTH => linear_settings%north
148  this%RELAX => linear_settings%relax
149  this%LEVEL => linear_settings%level
150  this%DROPTOL => linear_settings%droptol
151  !
152  ! -- SET POINTERS TO SOLUTION STORAGE
153  this%IPRIMS => iprims
154  this%NEQ => neq
155  call matrix%get_aij(this%IA, this%JA, this%AMAT)
156  call mem_allocate(this%NJA, 'NJA', this%memoryPath)
157  this%NJA = size(this%AMAT)
158  this%RHS => rhs
159  this%X => x
160  !
161  ! -- ALLOCATE SCALAR VARIABLES
162  call this%allocate_scalars()
163  !
164  ! -- initialize iout
165  this%iout = iout
166  !
167  ! -- DEFAULT VALUES
168  this%IPC = 0
169  !
170  this%IACPC = 0
171  !
172  ! -- PRINT A MESSAGE IDENTIFYING IMSLINEAR SOLVER PACKAGE
173  write (iout, 2000)
174 02000 FORMAT(1x, /1x, 'IMSLINEAR -- UNSTRUCTURED LINEAR SOLUTION', &
175  ' PACKAGE, VERSION 8, 04/28/2017')
176  !
177  ! -- DETERMINE PRECONDITIONER
178  IF (this%LEVEL > 0 .OR. this%DROPTOL > dzero) THEN
179  this%IPC = 3
180  ELSE
181  this%IPC = 1
182  END IF
183  IF (this%RELAX > dzero) THEN
184  this%IPC = this%IPC + 1
185  END IF
186  !
187  ! -- ERROR CHECKING FOR OPTIONS
188  IF (this%ISCL < 0) this%ISCL = 0
189  IF (this%ISCL > 2) THEN
190  WRITE (errmsg, '(A)') 'IMSLINEAR7AR ISCL MUST BE <= 2'
191  call store_error(errmsg)
192  END IF
193  IF (this%IORD < 0) this%IORD = 0
194  IF (this%IORD > 2) THEN
195  WRITE (errmsg, '(A)') 'IMSLINEAR7AR IORD MUST BE <= 2'
196  call store_error(errmsg)
197  END IF
198  IF (this%NORTH < 0) THEN
199  WRITE (errmsg, '(A)') 'IMSLINEAR7AR NORTH MUST >= 0'
200  call store_error(errmsg)
201  END IF
202  IF (this%RCLOSE == dzero) THEN
203  IF (this%ICNVGOPT /= 3) THEN
204  WRITE (errmsg, '(A)') 'IMSLINEAR7AR RCLOSE MUST > 0.0'
205  call store_error(errmsg)
206  END IF
207  END IF
208  IF (this%RELAX < dzero) THEN
209  WRITE (errmsg, '(A)') 'IMSLINEAR7AR RELAX MUST BE >= 0.0'
210  call store_error(errmsg)
211  END IF
212  IF (this%RELAX > done) THEN
213  WRITE (errmsg, '(A)') 'IMSLINEAR7AR RELAX MUST BE <= 1.0'
214  call store_error(errmsg)
215  END IF
216  !
217  ! -- INITIALIZE IMSLINEAR VARIABLES
218  this%NITERC = 0
219  !
220  ! -- ALLOCATE AND INITIALIZE MEMORY FOR IMSLINEAR
221  iscllen = 1
222  IF (this%ISCL .NE. 0) iscllen = neq
223  CALL mem_allocate(this%DSCALE, iscllen, 'DSCALE', trim(this%memoryPath))
224  CALL mem_allocate(this%DSCALE2, iscllen, 'DSCALE2', trim(this%memoryPath))
225  !
226  ! -- determine dimensions for preconditing arrays
227  call ims_calc_pcdims(this%NEQ, this%NJA, this%IA, this%LEVEL, this%IPC, &
228  this%NIAPC, this%NJAPC, this%NJLU, this%NJW, this%NWLU)
229  !
230  ! -- ALLOCATE BASE PRECONDITIONER VECTORS
231  CALL mem_allocate(this%IAPC, this%NIAPC + 1, 'IAPC', trim(this%memoryPath))
232  CALL mem_allocate(this%JAPC, this%NJAPC, 'JAPC', trim(this%memoryPath))
233  CALL mem_allocate(this%APC, this%NJAPC, 'APC', trim(this%memoryPath))
234  !
235  ! -- ALLOCATE MEMORY FOR ILU0 AND MILU0 NON-ZERO ROW ENTRY VECTOR
236  CALL mem_allocate(this%IW, this%NIAPC, 'IW', trim(this%memoryPath))
237  CALL mem_allocate(this%W, this%NIAPC, 'W', trim(this%memoryPath))
238  !
239  ! -- ALLOCATE MEMORY FOR ILUT VECTORS
240  CALL mem_allocate(this%JLU, this%NJLU, 'JLU', trim(this%memoryPath))
241  CALL mem_allocate(this%JW, this%NJW, 'JW', trim(this%memoryPath))
242  CALL mem_allocate(this%WLU, this%NWLU, 'WLU', trim(this%memoryPath))
243  !
244  ! -- GENERATE IAPC AND JAPC FOR ILU0 AND MILU0
245  IF (this%IPC == 1 .OR. this%IPC == 2) THEN
246  CALL ims_base_pccrs(this%NEQ, this%NJA, this%IA, this%JA, &
247  this%IAPC, this%JAPC)
248  END IF
249  !
250  ! -- ALLOCATE SPACE FOR PERMUTATION VECTOR
251  i0 = 1
252  iolen = 1
253  IF (this%IORD .NE. 0) THEN
254  i0 = this%NEQ
255  iolen = this%NJA
256  END IF
257  CALL mem_allocate(this%LORDER, i0, 'LORDER', trim(this%memoryPath))
258  CALL mem_allocate(this%IORDER, i0, 'IORDER', trim(this%memoryPath))
259  CALL mem_allocate(this%IARO, i0 + 1, 'IARO', trim(this%memoryPath))
260  CALL mem_allocate(this%JARO, iolen, 'JARO', trim(this%memoryPath))
261  CALL mem_allocate(this%ARO, iolen, 'ARO', trim(this%memoryPath))
262  !
263  ! -- ALLOCATE WORKING VECTORS FOR IMSLINEAR SOLVER
264  CALL mem_allocate(this%ID, this%NEQ, 'ID', trim(this%memoryPath))
265  CALL mem_allocate(this%D, this%NEQ, 'D', trim(this%memoryPath))
266  CALL mem_allocate(this%P, this%NEQ, 'P', trim(this%memoryPath))
267  CALL mem_allocate(this%Q, this%NEQ, 'Q', trim(this%memoryPath))
268  CALL mem_allocate(this%Z, this%NEQ, 'Z', trim(this%memoryPath))
269  !
270  ! -- ALLOCATE MEMORY FOR BCGS WORKING ARRAYS
271  this%NIABCGS = 1
272  IF (this%ILINMETH == 2) THEN
273  this%NIABCGS = this%NEQ
274  END IF
275  CALL mem_allocate(this%T, this%NIABCGS, 'T', trim(this%memoryPath))
276  CALL mem_allocate(this%V, this%NIABCGS, 'V', trim(this%memoryPath))
277  CALL mem_allocate(this%DHAT, this%NIABCGS, 'DHAT', trim(this%memoryPath))
278  CALL mem_allocate(this%PHAT, this%NIABCGS, 'PHAT', trim(this%memoryPath))
279  CALL mem_allocate(this%QHAT, this%NIABCGS, 'QHAT', trim(this%memoryPath))
280  !
281  ! -- INITIALIZE IMSLINEAR VECTORS
282  DO n = 1, iscllen
283  this%DSCALE(n) = done
284  this%DSCALE2(n) = done
285  END DO
286  DO n = 1, this%NJAPC
287  this%APC(n) = dzero
288  END DO
289  !
290  ! -- WORKING VECTORS
291  DO n = 1, this%NEQ
292  this%ID(n) = izero
293  this%D(n) = dzero
294  this%P(n) = dzero
295  this%Q(n) = dzero
296  this%Z(n) = dzero
297  END DO
298  DO n = 1, this%NIAPC
299  this%IW(n) = izero
300  this%W(n) = dzero
301  END DO
302  !
303  ! -- BCGS WORKING VECTORS
304  DO n = 1, this%NIABCGS
305  this%T(n) = dzero
306  this%V(n) = dzero
307  this%DHAT(n) = dzero
308  this%PHAT(n) = dzero
309  this%QHAT(n) = dzero
310  END DO
311  !
312  ! -- ILUT AND MILUT WORKING VECTORS
313  DO n = 1, this%NJLU
314  this%JLU(n) = izero
315  END DO
316  DO n = 1, this%NJW
317  this%JW(n) = izero
318  END DO
319  DO n = 1, this%NWLU
320  this%WLU(n) = dzero
321  END DO
322  !
323  ! -- REORDERING VECTORS
324  DO n = 1, i0 + 1
325  this%IARO(n) = izero
326  END DO
327  DO n = 1, iolen
328  this%JARO(n) = izero
329  this%ARO(n) = dzero
330  END DO
331  !
332  ! -- REVERSE CUTHILL MCKEE AND MINIMUM DEGREE ORDERING
333  IF (this%IORD .NE. 0) THEN
334  CALL ims_base_calc_order(this%IORD, this%NEQ, this%NJA, this%IA, &
335  this%JA, this%LORDER, this%IORDER)
336  END IF
337  !
338  ! -- ALLOCATE MEMORY FOR STORING ITERATION CONVERGENCE DATA
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
subroutine, public store_error(msg, terminate)
Store an error message.
Definition: Sim.f90:92
integer(i4b) function, public count_errors()
Return number of errors.
Definition: Sim.f90:59
subroutine, public deprecation_warning(cblock, cvar, cver, endmsg, iunit)
Store deprecation warning message.
Definition: Sim.f90:256
Here is the call graph for this function:

◆ imslinear_da()

subroutine imslinearmodule::imslinear_da ( class(imslineardatatype), intent(inout)  this)

Deallocate linear accelerator memory.

Parameters
[in,out]thislinear datatype instance

Definition at line 495 of file ImsLinear.f90.

496  ! -- modules
498  ! -- dummy variables
499  class(ImsLinearDataType), intent(inout) :: this !< linear datatype instance
500  !
501  ! -- arrays
502  call mem_deallocate(this%dscale)
503  call mem_deallocate(this%dscale2)
504  call mem_deallocate(this%iapc)
505  call mem_deallocate(this%japc)
506  call mem_deallocate(this%apc)
507  call mem_deallocate(this%iw)
508  call mem_deallocate(this%w)
509  call mem_deallocate(this%jlu)
510  call mem_deallocate(this%jw)
511  call mem_deallocate(this%wlu)
512  call mem_deallocate(this%lorder)
513  call mem_deallocate(this%iorder)
514  call mem_deallocate(this%iaro)
515  call mem_deallocate(this%jaro)
516  call mem_deallocate(this%aro)
517  call mem_deallocate(this%id)
518  call mem_deallocate(this%d)
519  call mem_deallocate(this%p)
520  call mem_deallocate(this%q)
521  call mem_deallocate(this%z)
522  call mem_deallocate(this%t)
523  call mem_deallocate(this%v)
524  call mem_deallocate(this%dhat)
525  call mem_deallocate(this%phat)
526  call mem_deallocate(this%qhat)
527  !
528  ! -- scalars
529  call mem_deallocate(this%iout)
530  call mem_deallocate(this%ipc)
531  call mem_deallocate(this%iacpc)
532  call mem_deallocate(this%niterc)
533  call mem_deallocate(this%niabcgs)
534  call mem_deallocate(this%niapc)
535  call mem_deallocate(this%njapc)
536  call mem_deallocate(this%epfact)
537  call mem_deallocate(this%l2norm0)
538  call mem_deallocate(this%njlu)
539  call mem_deallocate(this%njw)
540  call mem_deallocate(this%nwlu)
541  call mem_deallocate(this%NJA)
542  !
543  ! -- nullify pointers
544  nullify (this%iprims)
545  nullify (this%neq)
546  nullify (this%nja)
547  nullify (this%ia)
548  nullify (this%ja)
549  nullify (this%amat)
550  nullify (this%rhs)
551  nullify (this%x)

◆ imslinear_set_input()

subroutine imslinearmodule::imslinear_set_input ( class(imslineardatatype), intent(inout)  this,
integer(i4b), intent(in)  IFDPARAM 
)

Set default linear accelerator settings.

Parameters
[in,out]thisImsLinearDataType instance
[in]ifdparamcomplexity option

Definition at line 559 of file ImsLinear.f90.

560  ! -- dummy variables
561  CLASS(ImsLinearDataType), INTENT(INOUT) :: this !< ImsLinearDataType instance
562  integer(I4B), INTENT(IN) :: IFDPARAM !< complexity option
563  ! -- code
564  SELECT CASE (ifdparam)
565  !
566  ! -- Simple option
567  CASE (1)
568  this%ITER1 = 50
569  this%ILINMETH = 1
570  this%IPC = 1
571  this%ISCL = 0
572  this%IORD = 0
573  this%DVCLOSE = dem3
574  this%RCLOSE = dem1
575  this%RELAX = dzero
576  this%LEVEL = 0
577  this%DROPTOL = dzero
578  this%NORTH = 0
579  !
580  ! -- Moderate
581  CASE (2)
582  this%ITER1 = 100
583  this%ILINMETH = 2
584  this%IPC = 2
585  this%ISCL = 0
586  this%IORD = 0
587  this%DVCLOSE = dem2
588  this%RCLOSE = dem1
589  this%RELAX = 0.97d0
590  this%LEVEL = 0
591  this%DROPTOL = dzero
592  this%NORTH = 0
593  !
594  ! -- Complex
595  CASE (3)
596  this%ITER1 = 500
597  this%ILINMETH = 2
598  this%IPC = 3
599  this%ISCL = 0
600  this%IORD = 0
601  this%DVCLOSE = dem1
602  this%RCLOSE = dem1
603  this%RELAX = dzero
604  this%LEVEL = 5
605  this%DROPTOL = dem4
606  this%NORTH = 2
607  END SELECT

◆ imslinear_summary()

subroutine imslinearmodule::imslinear_summary ( class(imslineardatatype), intent(inout)  this,
integer(i4b), intent(in)  mxiter 
)

Write summary of linear accelerator settings.

Parameters
[in,out]thisImsLinearDataType instance
[in]mxitermaximum number of outer iterations

Definition at line 346 of file ImsLinear.f90.

347  ! -- dummy variables
348  class(ImsLinearDataType), intent(inout) :: this !< ImsLinearDataType instance
349  integer(I4B), intent(in) :: mxiter !< maximum number of outer iterations
350  ! -- local variables
351  CHARACTER(LEN=10) :: clin(0:2)
352  CHARACTER(LEN=31) :: clintit(0:2)
353  CHARACTER(LEN=20) :: cipc(0:4)
354  CHARACTER(LEN=20) :: cscale(0:2)
355  CHARACTER(LEN=25) :: corder(0:2)
356  CHARACTER(LEN=16), DIMENSION(0:4) :: ccnvgopt
357  CHARACTER(LEN=15) :: clevel
358  CHARACTER(LEN=15) :: cdroptol
359  integer(I4B) :: i
360  integer(I4B) :: j
361  ! -- data
362  DATA clin/'UNKNOWN ', &
363  &'CG ', &
364  &'BCGS '/
365  DATA clintit/' UNKNOWN ', &
366  &' CONJUGATE-GRADIENT ', &
367  &'BICONJUGATE-GRADIENT STABILIZED'/
368  DATA cipc/'UNKNOWN ', &
369  &'INCOMPLETE LU ', &
370  &'MOD. INCOMPLETE LU ', &
371  &'INCOMPLETE LUT ', &
372  &'MOD. INCOMPLETE LUT '/
373  DATA cscale/'NO SCALING ', &
374  &'SYMMETRIC SCALING ', &
375  &'L2 NORM SCALING '/
376  DATA corder/'ORIGINAL ORDERING ', &
377  &'RCM ORDERING ', &
378  &'MINIMUM DEGREE ORDERING '/
379  DATA ccnvgopt/'INFINITY NORM ', &
380  &'INFINITY NORM S ', &
381  &'L2 NORM ', &
382  &'RELATIVE L2NORM ', &
383  &'L2 NORM W. REL. '/
384  ! -- formats
385 02010 FORMAT(1x, /, 7x, 'SOLUTION BY THE', 1x, a31, 1x, 'METHOD', &
386  /, 1x, 66('-'), /, &
387  ' MAXIMUM OF ', i0, ' CALLS OF SOLUTION ROUTINE', /, &
388  ' MAXIMUM OF ', i0, &
389  ' INTERNAL ITERATIONS PER CALL TO SOLUTION ROUTINE', /, &
390  ' LINEAR ACCELERATION METHOD =', 1x, a, /, &
391  ' MATRIX PRECONDITIONING TYPE =', 1x, a, /, &
392  ' MATRIX SCALING APPROACH =', 1x, a, /, &
393  ' MATRIX REORDERING APPROACH =', 1x, a, /, &
394  ' NUMBER OF ORTHOGONALIZATIONS =', 1x, i0, /, &
395  ' HEAD CHANGE CRITERION FOR CLOSURE =', e15.5, /, &
396  ' RESIDUAL CHANGE CRITERION FOR CLOSURE =', e15.5, /, &
397  ' RESIDUAL CONVERGENCE OPTION =', 1x, i0, /, &
398  ' RESIDUAL CONVERGENCE NORM =', 1x, a, /, &
399  ' RELAXATION FACTOR =', e15.5)
400 02015 FORMAT(' NUMBER OF LEVELS =', a15, /, &
401  ' DROP TOLERANCE =', a15, //)
402 2030 FORMAT(1x, a20, 1x, 6(i6, 1x))
403 2040 FORMAT(1x, 20('-'), 1x, 6(6('-'), 1x))
404 2050 FORMAT(1x, 62('-'),/) !
405 ! -- -----------------------------------------------------------
406  !
407  ! -- initialize clevel and cdroptol
408  clevel = ''
409  cdroptol = ''
410  !
411  ! -- write common variables to all linear accelerators
412  write (this%iout, 2010) &
413  clintit(this%ILINMETH), mxiter, this%ITER1, &
414  clin(this%ILINMETH), cipc(this%IPC), &
415  cscale(this%ISCL), corder(this%IORD), &
416  this%NORTH, this%DVCLOSE, this%RCLOSE, &
417  this%ICNVGOPT, ccnvgopt(this%ICNVGOPT), &
418  this%RELAX
419  if (this%level > 0) then
420  write (clevel, '(i15)') this%level
421  end if
422  if (this%droptol > dzero) then
423  write (cdroptol, '(e15.5)') this%droptol
424  end if
425  IF (this%level > 0 .or. this%droptol > dzero) THEN
426  write (this%iout, 2015) trim(adjustl(clevel)), &
427  trim(adjustl(cdroptol))
428  ELSE
429  write (this%iout, '(//)')
430  END IF
431 
432  if (this%iord /= 0) then
433  !
434  ! -- WRITE SUMMARY OF REORDERING INFORMATION TO LIST FILE
435  if (this%iprims == 2) then
436  DO i = 1, this%neq, 6
437  write (this%iout, 2030) 'ORIGINAL NODE :', &
438  (j, j=i, min(i + 5, this%neq))
439  write (this%iout, 2040)
440  write (this%iout, 2030) 'REORDERED INDEX :', &
441  (this%lorder(j), j=i, min(i + 5, this%neq))
442  write (this%iout, 2030) 'REORDERED NODE :', &
443  (this%iorder(j), j=i, min(i + 5, this%neq))
444  write (this%iout, 2050)
445  END DO
446  END IF
447  end if