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

This module contains the IMS linear accelerator subroutines. More...

Functions/Subroutines

subroutine ims_base_cg (ICNVG, ITMAX, INNERIT, NEQ, NJA, NIAPC, NJAPC, IPC, ICNVGOPT, NORTH, DVCLOSE, RCLOSE, L2NORM0, EPFACT, IA0, JA0, A0, IAPC, JAPC, APC, X, B, D, P, Q, Z, NJLU, IW, JLU, NCONV, CONVNMOD, CONVMODSTART, CACCEL, summary)
 @ brief Preconditioned Conjugate Gradient linear accelerator More...
 
subroutine ims_base_bcgs (ICNVG, ITMAX, INNERIT, NEQ, NJA, NIAPC, NJAPC, IPC, ICNVGOPT, NORTH, ISCL, DSCALE, DVCLOSE, RCLOSE, L2NORM0, EPFACT, IA0, JA0, A0, IAPC, JAPC, APC, X, B, D, P, Q, T, V, DHAT, PHAT, QHAT, NJLU, IW, JLU, NCONV, CONVNMOD, CONVMODSTART, CACCEL, summary)
 @ brief Preconditioned BiConjugate Gradient Stabilized linear accelerator More...
 
subroutine ims_base_calc_order (IORD, NEQ, NJA, IA, JA, LORDER, IORDER)
 @ brief Calculate LORDER AND IORDER More...
 
subroutine ims_base_scale (IOPT, ISCL, NEQ, NJA, IA, JA, AMAT, X, B, DSCALE, DSCALE2)
 @ brief Scale the coefficient matrix More...
 
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 More...
 
subroutine ims_base_pcjac (NJA, NEQ, AMAT, APC, IA, JA)
 @ brief Jacobi preconditioner More...
 
subroutine ims_base_jaca (NEQ, A, D1, D2)
 @ brief Apply the Jacobi preconditioner More...
 
subroutine ims_base_pcilu0 (NJA, NEQ, AMAT, IA, JA, APC, IAPC, JAPC, IW, W, RELAX, IPCFLAG, DELTA)
 @ brief Update the ILU0 preconditioner More...
 
subroutine ims_base_ilu0a (NJA, NEQ, APC, IAPC, JAPC, R, D)
 @ brief Apply the ILU0 and MILU0 preconditioners More...
 
subroutine ims_base_testcnvg (Icnvgopt, Icnvg, Iiter, Dvmax, Rmax, Rmax0, Epfact, Dvclose, Rclose)
 @ brief Test for solver convergence More...
 
subroutine ims_calc_pcdims (neq, nja, ia, level, ipc, niapc, njapc, njlu, njw, nwlu)
 
subroutine ims_base_pccrs (NEQ, NJA, IA, JA, IAPC, JAPC)
 @ brief Generate CRS pointers for the preconditioner More...
 
subroutine ims_base_isort (NVAL, IARRAY)
 In-place sorting for an integer array. More...
 
subroutine ims_base_residual (NEQ, NJA, X, B, D, A, IA, JA)
 Calculate residual. More...
 
real(dp) function ims_base_epfact (icnvgopt, kstp)
 Function returning EPFACT. More...
 

Variables

type(blockparsertype), private parser
 

Detailed Description

This module contains the IMS linear accelerator subroutines used by a MODFLOW 6 solution.

Function/Subroutine Documentation

◆ ims_base_bcgs()

subroutine imslinearbasemodule::ims_base_bcgs ( integer(i4b), intent(inout)  ICNVG,
integer(i4b), intent(in)  ITMAX,
integer(i4b), intent(inout)  INNERIT,
integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NIAPC,
integer(i4b), intent(in)  NJAPC,
integer(i4b), intent(in)  IPC,
integer(i4b), intent(in)  ICNVGOPT,
integer(i4b), intent(in)  NORTH,
integer(i4b), intent(in)  ISCL,
real(dp), dimension(neq), intent(in)  DSCALE,
real(dp), intent(in)  DVCLOSE,
real(dp), intent(in)  RCLOSE,
real(dp), intent(in)  L2NORM0,
real(dp), intent(in)  EPFACT,
integer(i4b), dimension(neq + 1), intent(in)  IA0,
integer(i4b), dimension(nja), intent(in)  JA0,
real(dp), dimension(nja), intent(in)  A0,
integer(i4b), dimension(niapc + 1), intent(in)  IAPC,
integer(i4b), dimension(njapc), intent(in)  JAPC,
real(dp), dimension(njapc), intent(in)  APC,
real(dp), dimension(neq), intent(inout)  X,
real(dp), dimension(neq), intent(in)  B,
real(dp), dimension(neq), intent(inout)  D,
real(dp), dimension(neq), intent(inout)  P,
real(dp), dimension(neq), intent(inout)  Q,
real(dp), dimension(neq), intent(inout)  T,
real(dp), dimension(neq), intent(inout)  V,
real(dp), dimension(neq), intent(inout)  DHAT,
real(dp), dimension(neq), intent(inout)  PHAT,
real(dp), dimension(neq), intent(inout)  QHAT,
integer(i4b), intent(in)  NJLU,
integer(i4b), dimension(niapc), intent(in)  IW,
integer(i4b), dimension(njlu), intent(in)  JLU,
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 
)

Apply the Preconditioned BiConjugate Gradient Stabilized linear accelerator to the current coefficient matrix, right-hand side, using the currentdependent-variable.

Parameters
[in,out]icnvgconvergence flag (1) non-convergence (0)
[in]itmaxmaximum number of inner iterations
[in,out]inneritinner iteration count
[in]neqnumber of equations
[in]njanumber of non-zero entries
[in]niapcpreconditioner number of rows
[in]njapcpreconditioner number of non-zero entries
[in]ipcpreconditioner option
[in]icnvgoptflow convergence criteria option
[in]northorthogonalization frequency
[in]isclscaling option
[in]dscalescaling vector
[in]dvclosedependent-variable closure criteria
[in]rcloseflow closure criteria
[in]l2norm0initial L-2 norm for system of equations
[in]epfactfactor for decreasing flow convergence criteria for subsequent Picard iterations
[in]ia0CRS row pointers
[in]ja0CRS column pointers
[in]a0coefficient matrix
[in]iapcpreconditioner CRS row pointers
[in]japcpreconditioner CRS column pointers
[in]apcpreconditioner matrix
[in,out]xdependent-variable vector
[in]bright-hand side vector
[in,out]dpreconditioner working vector
[in,out]ppreconditioner working vector
[in,out]qpreconditioner working vector
[in,out]tpreconditioner working vector
[in,out]vpreconditioner working vector
[in,out]dhatBCGS preconditioner working vector
[in,out]phatBCGS preconditioner working vector
[in,out]qhatBCGS preconditioner working vector
[in]njlupreconditioner length of JLU vector
[in]iwpreconditioner integer working vector
[in]jlupreconditioner JLU working vector
[in]nconvmaximum number of inner iterations in a time step (maxiter * maxinner)
[in]convnmodnumber of models in the solution
[in,out]convmodstartpointer to the start of each model in the convmod* arrays
[in,out]caccelconvergence string
[in]summaryConvergence summary report

Definition at line 249 of file ImsLinearBase.f90.

259  ! -- dummy variables
260  integer(I4B), INTENT(INOUT) :: ICNVG !< convergence flag (1) non-convergence (0)
261  integer(I4B), INTENT(IN) :: ITMAX !< maximum number of inner iterations
262  integer(I4B), INTENT(INOUT) :: INNERIT !< inner iteration count
263  integer(I4B), INTENT(IN) :: NEQ !< number of equations
264  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
265  integer(I4B), INTENT(IN) :: NIAPC !< preconditioner number of rows
266  integer(I4B), INTENT(IN) :: NJAPC !< preconditioner number of non-zero entries
267  integer(I4B), INTENT(IN) :: IPC !< preconditioner option
268  integer(I4B), INTENT(IN) :: ICNVGOPT !< flow convergence criteria option
269  integer(I4B), INTENT(IN) :: NORTH !< orthogonalization frequency
270  integer(I4B), INTENT(IN) :: ISCL !< scaling option
271  real(DP), DIMENSION(NEQ), INTENT(IN) :: DSCALE !< scaling vector
272  real(DP), INTENT(IN) :: DVCLOSE !< dependent-variable closure criteria
273  real(DP), INTENT(IN) :: RCLOSE !< flow closure criteria
274  real(DP), INTENT(IN) :: L2NORM0 !< initial L-2 norm for system of equations
275  real(DP), INTENT(IN) :: EPFACT !< factor for decreasing flow convergence criteria for subsequent Picard iterations
276  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA0 !< CRS row pointers
277  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA0 !< CRS column pointers
278  real(DP), DIMENSION(NJA), INTENT(IN) :: A0 !< coefficient matrix
279  integer(I4B), DIMENSION(NIAPC + 1), INTENT(IN) :: IAPC !< preconditioner CRS row pointers
280  integer(I4B), DIMENSION(NJAPC), INTENT(IN) :: JAPC !< preconditioner CRS column pointers
281  real(DP), DIMENSION(NJAPC), INTENT(IN) :: APC !< preconditioner matrix
282  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: X !< dependent-variable vector
283  real(DP), DIMENSION(NEQ), INTENT(IN) :: B !< right-hand side vector
284  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: D !< preconditioner working vector
285  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: P !< preconditioner working vector
286  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: Q !< preconditioner working vector
287  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: T !< preconditioner working vector
288  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: V !< preconditioner working vector
289  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: DHAT !< BCGS preconditioner working vector
290  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: PHAT !< BCGS preconditioner working vector
291  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: QHAT !< BCGS preconditioner working vector
292  ! -- ILUT dummy variables
293  integer(I4B), INTENT(IN) :: NJLU !< preconditioner length of JLU vector
294  integer(I4B), DIMENSION(NIAPC), INTENT(IN) :: IW !< preconditioner integer working vector
295  integer(I4B), DIMENSION(NJLU), INTENT(IN) :: JLU !< preconditioner JLU working vector
296  ! -- convergence information dummy variables
297  integer(I4B), INTENT(IN) :: NCONV !< maximum number of inner iterations in a time step (maxiter * maxinner)
298  integer(I4B), INTENT(IN) :: CONVNMOD !< number of models in the solution
299  integer(I4B), DIMENSION(CONVNMOD + 1), INTENT(INOUT) :: CONVMODSTART !< pointer to the start of each model in the convmod* arrays
300  character(len=31), DIMENSION(NCONV), INTENT(INOUT) :: CACCEL !< convergence string
301  type(ConvergenceSummaryType), pointer, intent(in) :: summary !< Convergence summary report
302  ! -- local variables
303  LOGICAL :: LORTH
304  logical :: lsame
305  character(len=15) :: cval1, cval2
306  integer(I4B) :: n
307  integer(I4B) :: iiter
308  integer(I4B) :: xloc, rloc
309  integer(I4B) :: im, im0, im1
310  real(DP) :: ddot
311  real(DP) :: tv
312  real(DP) :: deltax
313  real(DP) :: rmax
314  real(DP) :: l2norm
315  real(DP) :: rcnvg
316  real(DP) :: alpha, alpha0
317  real(DP) :: beta
318  real(DP) :: rho, rho0
319  real(DP) :: omega, omega0
320  real(DP) :: numerator, denominator
321  !
322  ! -- initialize local variables
323  innerit = 0
324  alpha = dzero
325  alpha0 = dzero
326  beta = dzero
327  rho = dzero
328  rho0 = dzero
329  omega = dzero
330  omega0 = dzero
331  !
332  ! -- SAVE INITIAL RESIDUAL
333  DO n = 1, neq
334  dhat(n) = d(n)
335  END DO
336  !
337  ! -- INNER ITERATION
338  inner: DO iiter = 1, itmax
339  innerit = innerit + 1
340  summary%iter_cnt = summary%iter_cnt + 1
341  !
342  ! -- CALCULATE rho
343  rho = ddot(neq, dhat, 1, d, 1)
344  !
345  ! -- COMPUTE DIRECTIONAL VECTORS
346  IF (iiter == 1) THEN
347  DO n = 1, neq
348  p(n) = d(n)
349  END DO
350  ELSE
351  beta = (rho / rho0) * (alpha0 / omega0)
352  DO n = 1, neq
353  p(n) = d(n) + beta * (p(n) - omega0 * v(n))
354  END DO
355  END IF
356  !
357  ! -- APPLY PRECONDITIONER TO UPDATE PHAT
358  SELECT CASE (ipc)
359  !
360  ! -- ILU0 AND MILU0
361  CASE (1, 2)
362  CALL ims_base_ilu0a(nja, neq, apc, iapc, japc, p, phat)
363  !
364  ! -- ILUT AND MILUT
365  CASE (3, 4)
366  CALL lusol(neq, p, phat, apc, jlu, iw)
367  END SELECT
368  !
369  ! -- COMPUTE ITERATES
370  !
371  ! -- UPDATE V WITH A AND PHAT
372  call amux(neq, phat, v, a0, ja0, ia0)
373  !
374  ! -- UPDATE alpha WITH DHAT AND V
375  denominator = ddot(neq, dhat, 1, v, 1)
376  denominator = denominator + sign(dprec, denominator)
377  alpha = rho / denominator
378  !
379  ! -- UPDATE Q
380  DO n = 1, neq
381  q(n) = d(n) - alpha * v(n)
382  END DO
383  !
384  ! ! -- CALCULATE INFINITY NORM OF Q - TEST FOR TERMINATION
385  ! ! TERMINATE IF rmax IS LESS THAN MACHINE PRECISION (DPREC)
386  ! rmax = DZERO
387  ! DO n = 1, NEQ
388  ! tv = Q(n)
389  ! IF (ISCL.NE.0 ) tv = tv / DSCALE(n)
390  ! IF (ABS(tv) > ABS(rmax) ) rmax = tv
391  ! END DO
392  ! IF (ABS(rmax).LE.DPREC) THEN
393  ! deltax = DZERO
394  ! DO n = 1, NEQ
395  ! tv = alpha * PHAT(n)
396  ! IF (ISCL.NE.0) THEN
397  ! tv = tv * DSCALE(n)
398  ! END IF
399  ! X(n) = X(n) + tv
400  ! IF (ABS(tv) > ABS(deltax) ) deltax = tv
401  ! END DO
402  ! CALL IMSLINEARSUB_TESTCNVG(ICNVGOPT, ICNVG, INNERIT, &
403  ! deltax, rmax, &
404  ! rmax, EPFACT, DVCLOSE, RCLOSE )
405  ! IF (ICNVG.NE.0 ) EXIT INNER
406  ! END IF
407  !
408  ! -- APPLY PRECONDITIONER TO UPDATE QHAT
409  SELECT CASE (ipc)
410  !
411  ! -- ILU0 AND MILU0
412  CASE (1, 2)
413  CALL ims_base_ilu0a(nja, neq, apc, iapc, japc, q, qhat)
414  !
415  ! -- ILUT AND MILUT
416  CASE (3, 4)
417  CALL lusol(neq, q, qhat, apc, jlu, iw)
418  END SELECT
419  !
420  ! -- UPDATE T WITH A AND QHAT
421  call amux(neq, qhat, t, a0, ja0, ia0)
422  !
423  ! -- UPDATE omega
424  numerator = ddot(neq, t, 1, q, 1)
425  denominator = ddot(neq, t, 1, t, 1)
426  denominator = denominator + sign(dprec, denominator)
427  omega = numerator / denominator
428  !
429  ! -- UPDATE X AND RESIDUAL
430  deltax = dzero
431  rmax = dzero
432  l2norm = dzero
433  DO im = 1, convnmod
434  summary%dvmax(im) = dzero
435  summary%rmax(im) = dzero
436  END DO
437  im = 1
438  im0 = convmodstart(1)
439  im1 = convmodstart(2)
440  DO n = 1, neq
441  !
442  ! -- determine current model index
443  if (n == im1) then
444  im = im + 1
445  im0 = convmodstart(im)
446  im1 = convmodstart(im + 1)
447  end if
448  !
449  ! -- X AND DX
450  tv = alpha * phat(n) + omega * qhat(n)
451  x(n) = x(n) + tv
452  IF (iscl .NE. 0) THEN
453  tv = tv * dscale(n)
454  END IF
455  IF (abs(tv) > abs(deltax)) THEN
456  deltax = tv
457  xloc = n
458  END IF
459  IF (abs(tv) > abs(summary%dvmax(im))) THEN
460  summary%dvmax(im) = tv
461  summary%locdv(im) = n
462  END IF
463  !
464  ! -- RESIDUAL
465  tv = q(n) - omega * t(n)
466  d(n) = tv
467  IF (iscl .NE. 0) THEN
468  tv = tv / dscale(n)
469  END IF
470  IF (abs(tv) > abs(rmax)) THEN
471  rmax = tv
472  rloc = n
473  END IF
474  IF (abs(tv) > abs(summary%rmax(im))) THEN
475  summary%rmax(im) = tv
476  summary%locr(im) = n
477  END IF
478  l2norm = l2norm + tv * tv
479  END DO
480  l2norm = sqrt(l2norm)
481  !
482  ! -- SAVE SOLVER convergence information dummy variables
483  IF (nconv > 1) THEN !<
484  n = summary%iter_cnt
485  WRITE (cval1, '(g15.7)') alpha
486  WRITE (cval2, '(g15.7)') omega
487  caccel(n) = trim(adjustl(cval1))//','//trim(adjustl(cval2))
488  summary%itinner(n) = iiter
489  DO im = 1, convnmod
490  summary%convdvmax(im, n) = summary%dvmax(im)
491  summary%convlocdv(im, n) = summary%locdv(im)
492  summary%convrmax(im, n) = summary%rmax(im)
493  summary%convlocr(im, n) = summary%locr(im)
494  END DO
495  END IF
496  !
497  ! -- TEST FOR SOLVER CONVERGENCE
498  IF (icnvgopt == 2 .OR. icnvgopt == 3 .OR. icnvgopt == 4) THEN
499  rcnvg = l2norm
500  ELSE
501  rcnvg = rmax
502  END IF
503  CALL ims_base_testcnvg(icnvgopt, icnvg, innerit, &
504  deltax, rcnvg, &
505  l2norm0, epfact, dvclose, rclose)
506  !
507  ! -- CHECK FOR EXACT SOLUTION
508  IF (rcnvg == dzero) icnvg = 1
509  !
510  ! -- CHECK FOR STANDARD CONVERGENCE
511  IF (icnvg .NE. 0) EXIT inner
512  !
513  ! -- CHECK THAT CURRENT AND PREVIOUS rho, alpha, AND omega ARE
514  ! DIFFERENT
515  lsame = is_close(rho, rho0)
516  IF (lsame) THEN
517  EXIT inner
518  END IF
519  lsame = is_close(alpha, alpha0)
520  IF (lsame) THEN
521  EXIT inner
522  END IF
523  lsame = is_close(omega, omega0)
524  IF (lsame) THEN
525  EXIT inner
526  END IF
527  !
528  ! -- RECALCULATE THE RESIDUAL
529  IF (north > 0) THEN
530  lorth = mod(iiter + 1, north) == 0
531  IF (lorth) THEN
532  call ims_base_residual(neq, nja, x, b, d, a0, ia0, ja0)
533  END IF
534  END IF
535  !
536  ! -- exit inner if rho or omega are zero
537  if (rho * omega == dzero) then
538  exit inner
539  end if
540  !
541  ! -- SAVE CURRENT INNER ITERATES
542  rho0 = rho
543  alpha0 = alpha
544  omega0 = omega
545  END DO inner
546  !
547  ! -- RESET ICNVG
548  IF (icnvg < 0) icnvg = 0
real(kind=8) function ddot(n, dx, incx, dy, incy)
Definition: blas1_d.f90:296
subroutine lusol(n, y, x, alu, jlu, ju)
Definition: ilut.f90:466
subroutine amux(n, x, y, a, ja, ia)
Definition: sparsekit.f90:2
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_calc_order()

subroutine imslinearbasemodule::ims_base_calc_order ( integer(i4b), intent(in)  IORD,
integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA,
integer(i4b), dimension(neq), intent(inout)  LORDER,
integer(i4b), dimension(neq), intent(inout)  IORDER 
)

Calculate LORDER and IORDER for reordering.

Parameters
[in]iordreordering option
[in]neqnumber of rows
[in]njanumber of non-zero entries
[in]iarow pointer
[in]jacolumn pointer
[in,out]lorderreorder vector
[in,out]iorderinverse of reorder vector

Definition at line 556 of file ImsLinearBase.f90.

557  ! -- modules
559  ! -- dummy variables
560  integer(I4B), INTENT(IN) :: IORD !< reordering option
561  integer(I4B), INTENT(IN) :: NEQ !< number of rows
562  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
563  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< row pointer
564  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< column pointer
565  integer(I4B), DIMENSION(NEQ), INTENT(INOUT) :: LORDER !< reorder vector
566  integer(I4B), DIMENSION(NEQ), INTENT(INOUT) :: IORDER !< inverse of reorder vector
567  ! -- local variables
568  character(len=LINELENGTH) :: errmsg
569  integer(I4B) :: n
570  integer(I4B) :: nsp
571  integer(I4B), DIMENSION(:), ALLOCATABLE :: iwork0
572  integer(I4B), DIMENSION(:), ALLOCATABLE :: iwork1
573  integer(I4B) :: iflag
574  !
575  ! -- initialize lorder and iorder
576  DO n = 1, neq
577  lorder(n) = izero
578  iorder(n) = izero
579  END DO
580  ! ALLOCATE (iwork0(NEQ))
581  SELECT CASE (iord)
582  CASE (1)
583  CALL genrcm(neq, nja, ia, ja, lorder)
584  CASE (2)
585  nsp = 3 * neq + 4 * nja
586  allocate (iwork0(neq))
587  allocate (iwork1(nsp))
588  CALL ims_odrv(neq, nja, nsp, ia, ja, lorder, iwork0, &
589  iwork1, iflag)
590  IF (iflag .NE. 0) THEN
591  write (errmsg, '(A,1X,A)') &
592  'IMSLINEARSUB_CALC_ORDER error creating minimum degree ', &
593  'order permutation '
594  call store_error(errmsg)
595  END IF
596  !
597  ! -- DEALLOCATE TEMPORARY STORAGE
598  deallocate (iwork0, iwork1)
599  END SELECT
600  !
601  ! -- GENERATE INVERSE OF LORDER
602  DO n = 1, neq
603  iorder(lorder(n)) = n
604  END DO
605  !
606  ! -- terminate if errors occurred
607  if (count_errors() > 0) then
608  call parser%StoreErrorUnit()
609  end if
This module contains simulation methods.
Definition: Sim.f90:10
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 genrcm(node_num, adj_num, adj_row, adj, perm)
Definition: rcm.f90:1004
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_cg()

subroutine imslinearbasemodule::ims_base_cg ( integer(i4b), intent(inout)  ICNVG,
integer(i4b), intent(in)  ITMAX,
integer(i4b), intent(inout)  INNERIT,
integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NIAPC,
integer(i4b), intent(in)  NJAPC,
integer(i4b), intent(in)  IPC,
integer(i4b), intent(in)  ICNVGOPT,
integer(i4b), intent(in)  NORTH,
real(dp), intent(in)  DVCLOSE,
real(dp), intent(in)  RCLOSE,
real(dp), intent(in)  L2NORM0,
real(dp), intent(in)  EPFACT,
integer(i4b), dimension(neq + 1), intent(in)  IA0,
integer(i4b), dimension(nja), intent(in)  JA0,
real(dp), dimension(nja), intent(in)  A0,
integer(i4b), dimension(niapc + 1), intent(in)  IAPC,
integer(i4b), dimension(njapc), intent(in)  JAPC,
real(dp), dimension(njapc), intent(in)  APC,
real(dp), dimension(neq), intent(inout)  X,
real(dp), dimension(neq), intent(inout)  B,
real(dp), dimension(neq), intent(inout)  D,
real(dp), dimension(neq), intent(inout)  P,
real(dp), dimension(neq), intent(inout)  Q,
real(dp), dimension(neq), intent(inout)  Z,
integer(i4b), intent(in)  NJLU,
integer(i4b), dimension(niapc), intent(in)  IW,
integer(i4b), dimension(njlu), intent(in)  JLU,
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 
)

Apply the Preconditioned Conjugate Gradient linear accelerator to the current coefficient matrix, right-hand side, using the current dependent-variable.

Parameters
[in,out]icnvgconvergence flag (1) non-convergence (0)
[in]itmaxmaximum number of inner iterations
[in,out]inneritinner iteration count
[in]neqnumber of equations
[in]njanumber of non-zero entries
[in]niapcpreconditioner number of rows
[in]njapcpreconditioner number of non-zero entries
[in]ipcpreconditioner option
[in]icnvgoptflow convergence criteria option
[in]northorthogonalization frequency
[in]dvclosedependent-variable closure criteria
[in]rcloseflow closure criteria
[in]l2norm0initial L-2 norm for system of equations
[in]epfactfactor for decreasing flow convergence criteria for subsequent Picard iterations
[in]ia0CRS row pointers
[in]ja0CRS column pointers
[in]a0coefficient matrix
[in]iapcpreconditioner CRS row pointers
[in]japcpreconditioner CRS column pointers
[in]apcpreconditioner matrix
[in,out]xdependent-variable vector
[in,out]bright-hand side vector
[in,out]dworking vector
[in,out]pworking vector
[in,out]qworking vector
[in,out]zworking vector
[in]njlupreconditioner length of JLU vector
[in]iwpreconditioner integer working vector
[in]jlupreconditioner JLU working vector
[in]nconvmaximum number of inner iterations in a time step (maxiter * maxinner)
[in]convnmodnumber of models in the solution
[in,out]convmodstartpointer to the start of each model in the convmod* arrays
[in,out]caccelconvergence string
[in]summaryConvergence summary report

Definition at line 30 of file ImsLinearBase.f90.

39  ! -- dummy variables
40  integer(I4B), INTENT(INOUT) :: ICNVG !< convergence flag (1) non-convergence (0)
41  integer(I4B), INTENT(IN) :: ITMAX !< maximum number of inner iterations
42  integer(I4B), INTENT(INOUT) :: INNERIT !< inner iteration count
43  integer(I4B), INTENT(IN) :: NEQ !< number of equations
44  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
45  integer(I4B), INTENT(IN) :: NIAPC !< preconditioner number of rows
46  integer(I4B), INTENT(IN) :: NJAPC !< preconditioner number of non-zero entries
47  integer(I4B), INTENT(IN) :: IPC !< preconditioner option
48  integer(I4B), INTENT(IN) :: ICNVGOPT !< flow convergence criteria option
49  integer(I4B), INTENT(IN) :: NORTH !< orthogonalization frequency
50  real(DP), INTENT(IN) :: DVCLOSE !< dependent-variable closure criteria
51  real(DP), INTENT(IN) :: RCLOSE !< flow closure criteria
52  real(DP), INTENT(IN) :: L2NORM0 !< initial L-2 norm for system of equations
53  real(DP), INTENT(IN) :: EPFACT !< factor for decreasing flow convergence criteria for subsequent Picard iterations
54  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA0 !< CRS row pointers
55  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA0 !< CRS column pointers
56  real(DP), DIMENSION(NJA), INTENT(IN) :: A0 !< coefficient matrix
57  integer(I4B), DIMENSION(NIAPC + 1), INTENT(IN) :: IAPC !< preconditioner CRS row pointers
58  integer(I4B), DIMENSION(NJAPC), INTENT(IN) :: JAPC !< preconditioner CRS column pointers
59  real(DP), DIMENSION(NJAPC), INTENT(IN) :: APC !< preconditioner matrix
60  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: X !< dependent-variable vector
61  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: B !< right-hand side vector
62  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: D !< working vector
63  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: P !< working vector
64  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: Q !< working vector
65  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: Z !< working vector
66  ! -- ILUT dummy variables
67  integer(I4B), INTENT(IN) :: NJLU !< preconditioner length of JLU vector
68  integer(I4B), DIMENSION(NIAPC), INTENT(IN) :: IW !< preconditioner integer working vector
69  integer(I4B), DIMENSION(NJLU), INTENT(IN) :: JLU !< preconditioner JLU working vector
70  ! -- convergence information dummy variables dummy variables
71  integer(I4B), INTENT(IN) :: NCONV !< maximum number of inner iterations in a time step (maxiter * maxinner)
72  integer(I4B), INTENT(IN) :: CONVNMOD !< number of models in the solution
73  integer(I4B), DIMENSION(CONVNMOD + 1), INTENT(INOUT) :: CONVMODSTART !< pointer to the start of each model in the convmod* arrays
74  character(len=31), DIMENSION(NCONV), INTENT(INOUT) :: CACCEL !< convergence string
75  type(ConvergenceSummaryType), pointer, intent(in) :: summary !< Convergence summary report
76  ! -- local variables
77  LOGICAL :: lorth
78  logical :: lsame
79  character(len=31) :: cval
80  integer(I4B) :: n
81  integer(I4B) :: iiter
82  integer(I4B) :: xloc, rloc
83  integer(I4B) :: im, im0, im1
84  real(DP) :: ddot
85  real(DP) :: tv
86  real(DP) :: deltax
87  real(DP) :: rmax
88  real(DP) :: l2norm
89  real(DP) :: rcnvg
90  real(DP) :: denominator
91  real(DP) :: alpha, beta
92  real(DP) :: rho, rho0
93  !
94  ! -- initialize local variables
95  rho0 = dzero
96  rho = dzero
97  innerit = 0
98  !
99  ! -- INNER ITERATION
100  inner: DO iiter = 1, itmax
101  innerit = innerit + 1
102  summary%iter_cnt = summary%iter_cnt + 1
103  !
104  ! -- APPLY PRECONDITIONER
105  SELECT CASE (ipc)
106  !
107  ! -- ILU0 AND MILU0
108  CASE (1, 2)
109  CALL ims_base_ilu0a(nja, neq, apc, iapc, japc, d, z)
110  !
111  ! -- ILUT AND MILUT
112  CASE (3, 4)
113  CALL lusol(neq, d, z, apc, jlu, iw)
114  END SELECT
115  rho = ddot(neq, d, 1, z, 1)
116  !
117  ! -- COMPUTE DIRECTIONAL VECTORS
118  IF (iiter == 1) THEN
119  DO n = 1, neq
120  p(n) = z(n)
121  END DO
122  ELSE
123  beta = rho / rho0
124  DO n = 1, neq
125  p(n) = z(n) + beta * p(n)
126  END DO
127  END IF
128  !
129  ! -- COMPUTE ITERATES
130  !
131  ! -- UPDATE Q
132  call amux(neq, p, q, a0, ja0, ia0)
133  denominator = ddot(neq, p, 1, q, 1)
134  denominator = denominator + sign(dprec, denominator)
135  alpha = rho / denominator
136  !
137  ! -- UPDATE X AND RESIDUAL
138  deltax = dzero
139  rmax = dzero
140  l2norm = dzero
141  DO im = 1, convnmod
142  summary%locdv(im) = 0
143  summary%dvmax(im) = dzero
144  summary%locr(im) = 0
145  summary%rmax(im) = dzero
146  END DO
147  im = 1
148  im0 = convmodstart(1)
149  im1 = convmodstart(2)
150  DO n = 1, neq
151  !
152  ! -- determine current model index
153  if (n == im1) then
154  im = im + 1
155  im0 = convmodstart(im)
156  im1 = convmodstart(im + 1)
157  end if
158  !
159  ! -- identify deltax and rmax
160  tv = alpha * p(n)
161  x(n) = x(n) + tv
162  IF (abs(tv) > abs(deltax)) THEN
163  deltax = tv
164  xloc = n
165  END IF
166  IF (abs(tv) > abs(summary%dvmax(im))) THEN
167  summary%dvmax(im) = tv
168  summary%locdv(im) = n
169  END IF
170  tv = d(n)
171  tv = tv - alpha * q(n)
172  d(n) = tv
173  IF (abs(tv) > abs(rmax)) THEN
174  rmax = tv
175  rloc = n
176  END IF
177  IF (abs(tv) > abs(summary%rmax(im))) THEN
178  summary%rmax(im) = tv
179  summary%locr(im) = n
180  END IF
181  l2norm = l2norm + tv * tv
182  END DO
183  l2norm = sqrt(l2norm)
184  !
185  ! -- SAVE SOLVER convergence information dummy variables
186  IF (nconv > 1) THEN !<
187  n = summary%iter_cnt
188  WRITE (cval, '(g15.7)') alpha
189  caccel(n) = cval
190  summary%itinner(n) = iiter
191  DO im = 1, convnmod
192  summary%convlocdv(im, n) = summary%locdv(im)
193  summary%convlocr(im, n) = summary%locr(im)
194  summary%convdvmax(im, n) = summary%dvmax(im)
195  summary%convrmax(im, n) = summary%rmax(im)
196  END DO
197  END IF
198  !
199  ! -- TEST FOR SOLVER CONVERGENCE
200  IF (icnvgopt == 2 .OR. icnvgopt == 3 .OR. icnvgopt == 4) THEN
201  rcnvg = l2norm
202  ELSE
203  rcnvg = rmax
204  END IF
205  CALL ims_base_testcnvg(icnvgopt, icnvg, innerit, &
206  deltax, rcnvg, &
207  l2norm0, epfact, dvclose, rclose)
208  !
209  ! -- CHECK FOR EXACT SOLUTION
210  IF (rcnvg == dzero) icnvg = 1
211  !
212  ! -- CHECK FOR STANDARD CONVERGENCE
213  IF (icnvg .NE. 0) EXIT inner
214  !
215  ! -- CHECK THAT CURRENT AND PREVIOUS rho ARE DIFFERENT
216  lsame = is_close(rho, rho0)
217  IF (lsame) THEN
218  EXIT inner
219  END IF
220  !
221  ! -- RECALCULATE THE RESIDUAL
222  IF (north > 0) THEN
223  lorth = mod(iiter + 1, north) == 0
224  IF (lorth) THEN
225  call ims_base_residual(neq, nja, x, b, d, a0, ia0, ja0)
226  END IF
227  END IF
228  !
229  ! -- exit inner if rho is zero
230  if (rho == dzero) then
231  exit inner
232  end if
233  !
234  ! -- SAVE CURRENT INNER ITERATES
235  rho0 = rho
236  END DO inner
237  !
238  ! -- RESET ICNVG
239  IF (icnvg < 0) icnvg = 0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_epfact()

real(dp) function imslinearbasemodule::ims_base_epfact ( integer(i4b)  icnvgopt,
integer(i4b)  kstp 
)
Parameters
icnvgoptIMS convergence option
kstptime step number
Returns
factor for decreasing convergence criteria in subsequent Picard iterations

Definition at line 1316 of file ImsLinearBase.f90.

1317  integer(I4B) :: icnvgopt !< IMS convergence option
1318  integer(I4B) :: kstp !< time step number
1319  real(DP) :: epfact !< factor for decreasing convergence criteria in subsequent Picard iterations
1320 
1321  if (icnvgopt == 2) then
1322  if (kstp == 1) then
1323  epfact = 0.01
1324  else
1325  epfact = 0.10
1326  end if
1327  else if (icnvgopt == 4) then
1328  epfact = dem4
1329  else
1330  epfact = done
1331  end if
1332 
Here is the caller graph for this function:

◆ ims_base_ilu0a()

subroutine imslinearbasemodule::ims_base_ilu0a ( integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NEQ,
real(dp), dimension(nja), intent(in)  APC,
integer(i4b), dimension(neq + 1), intent(in)  IAPC,
integer(i4b), dimension(nja), intent(in)  JAPC,
real(dp), dimension(neq), intent(in)  R,
real(dp), dimension(neq), intent(inout)  D 
)

Apply the ILU0 and MILU0 preconditioners to the passed vector (R).

Parameters
[in]njanumber of non-zero entries
[in]neqnumber of equations
[in]apcILU0/MILU0 preconditioner matrix
[in]iapcILU0/MILU0 preconditioner CRS row pointers
[in]japcILU0/MILU0 preconditioner CRS column pointers
[in]rinput vector
[in,out]doutput vector after applying APC to R

Definition at line 1049 of file ImsLinearBase.f90.

1050  ! -- dummy variables
1051  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
1052  integer(I4B), INTENT(IN) :: NEQ !< number of equations
1053  real(DP), DIMENSION(NJA), INTENT(IN) :: APC !< ILU0/MILU0 preconditioner matrix
1054  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IAPC !< ILU0/MILU0 preconditioner CRS row pointers
1055  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JAPC !< ILU0/MILU0 preconditioner CRS column pointers
1056  real(DP), DIMENSION(NEQ), INTENT(IN) :: R !< input vector
1057  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: D !< output vector after applying APC to R
1058  ! -- local variables
1059  integer(I4B) :: ic0, ic1
1060  integer(I4B) :: iu
1061  integer(I4B) :: jcol
1062  integer(I4B) :: j, n
1063  real(DP) :: tv
1064  !
1065  ! -- FORWARD SOLVE - APC * D = R
1066  forward: DO n = 1, neq
1067  tv = r(n)
1068  ic0 = iapc(n)
1069  ic1 = iapc(n + 1) - 1
1070  iu = japc(n) - 1
1071  lower: DO j = ic0, iu
1072  jcol = japc(j)
1073  tv = tv - apc(j) * d(jcol)
1074  END DO lower
1075  d(n) = tv
1076  END DO forward
1077  !
1078  ! -- BACKWARD SOLVE - D = D / U
1079  backward: DO n = neq, 1, -1
1080  ic0 = iapc(n)
1081  ic1 = iapc(n + 1) - 1
1082  iu = japc(n)
1083  tv = d(n)
1084  upper: DO j = iu, ic1
1085  jcol = japc(j)
1086  tv = tv - apc(j) * d(jcol)
1087  END DO upper
1088  !
1089  ! -- COMPUTE D FOR DIAGONAL - D = D / U
1090  d(n) = tv * apc(n)
1091  END DO backward
Here is the caller graph for this function:

◆ ims_base_isort()

subroutine imslinearbasemodule::ims_base_isort ( integer(i4b), intent(in)  NVAL,
integer(i4b), dimension(nval), intent(inout)  IARRAY 
)

Subroutine sort an integer array in-place.

Parameters
[in]nvallength of the integer array
[in,out]iarrayinteger array to be sorted

Definition at line 1268 of file ImsLinearBase.f90.

1269  ! -- dummy variables
1270  integer(I4B), INTENT(IN) :: NVAL !< length of the integer array
1271  integer(I4B), DIMENSION(NVAL), INTENT(INOUT) :: IARRAY !< integer array to be sorted
1272  ! -- local variables
1273  integer(I4B) :: i, j, itemp
1274  ! -- code
1275  DO i = 1, nval - 1
1276  DO j = i + 1, nval
1277  if (iarray(i) > iarray(j)) then
1278  itemp = iarray(j)
1279  iarray(j) = iarray(i)
1280  iarray(i) = itemp
1281  END IF
1282  END DO
1283  END DO
Here is the caller graph for this function:

◆ ims_base_jaca()

subroutine imslinearbasemodule::ims_base_jaca ( integer(i4b), intent(in)  NEQ,
real(dp), dimension(neq), intent(in)  A,
real(dp), dimension(neq), intent(in)  D1,
real(dp), dimension(neq), intent(inout)  D2 
)

Apply the Jacobi preconditioner and return the resultant vector.

Parameters
[in]neqnumber of equations
[in]aJacobi preconditioner
[in]d1input vector
[in,out]d2resultant vector

Definition at line 907 of file ImsLinearBase.f90.

908  ! -- dummy variables
909  integer(I4B), INTENT(IN) :: NEQ !< number of equations
910  real(DP), DIMENSION(NEQ), INTENT(IN) :: A !< Jacobi preconditioner
911  real(DP), DIMENSION(NEQ), INTENT(IN) :: D1 !< input vector
912  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: D2 !< resultant vector
913  ! -- local variables
914  integer(I4B) :: n
915  real(DP) :: tv
916  ! -- code
917  DO n = 1, neq
918  tv = a(n) * d1(n)
919  d2(n) = tv
920  END DO

◆ ims_base_pccrs()

subroutine imslinearbasemodule::ims_base_pccrs ( integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA,
integer(i4b), dimension(neq + 1), intent(inout)  IAPC,
integer(i4b), dimension(nja), intent(inout)  JAPC 
)

Generate the CRS row and column pointers for the preconditioner. JAPC(1:NEQ) hHas the position of the upper entry for a row, JAPC(NEQ+1:NJA) is the column position for entry, APC(1:NEQ) is the preconditioned inverse of the diagonal, and APC(NEQ+1:NJA) are the preconditioned entries for off diagonals.

Definition at line 1207 of file ImsLinearBase.f90.

1209  ! -- dummy variables
1210  integer(I4B), INTENT(IN) :: NEQ !<
1211  integer(I4B), INTENT(IN) :: NJA !<
1212  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !<
1213  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !<
1214  integer(I4B), DIMENSION(NEQ + 1), INTENT(INOUT) :: IAPC !<
1215  integer(I4B), DIMENSION(NJA), INTENT(INOUT) :: JAPC !<
1216  ! -- local variables
1217  integer(I4B) :: n, j
1218  integer(I4B) :: i0, i1
1219  integer(I4B) :: nlen
1220  integer(I4B) :: ic, ip
1221  integer(I4B) :: jcol
1222  integer(I4B), DIMENSION(:), ALLOCATABLE :: iarr
1223  ! -- code
1224  ip = neq + 1
1225  DO n = 1, neq
1226  i0 = ia(n)
1227  i1 = ia(n + 1) - 1
1228  nlen = i1 - i0
1229  ALLOCATE (iarr(nlen))
1230  ic = 0
1231  DO j = i0, i1
1232  jcol = ja(j)
1233  IF (jcol == n) cycle
1234  ic = ic + 1
1235  iarr(ic) = jcol
1236  END DO
1237  CALL ims_base_isort(nlen, iarr)
1238  iapc(n) = ip
1239  DO j = 1, nlen
1240  jcol = iarr(j)
1241  japc(ip) = jcol
1242  ip = ip + 1
1243  END DO
1244  DEALLOCATE (iarr)
1245  END DO
1246  iapc(neq + 1) = nja + 1
1247  !
1248  ! -- POSITION OF THE FIRST UPPER ENTRY FOR ROW
1249  DO n = 1, neq
1250  i0 = iapc(n)
1251  i1 = iapc(n + 1) - 1
1252  japc(n) = iapc(n + 1)
1253  DO j = i0, i1
1254  jcol = japc(j)
1255  IF (jcol > n) THEN
1256  japc(n) = j
1257  EXIT
1258  END IF
1259  END DO
1260  END DO
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_pcilu0()

subroutine imslinearbasemodule::ims_base_pcilu0 ( integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NEQ,
real(dp), dimension(nja), intent(in)  AMAT,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA,
real(dp), dimension(nja), intent(inout)  APC,
integer(i4b), dimension(neq + 1), intent(inout)  IAPC,
integer(i4b), dimension(nja), intent(inout)  JAPC,
integer(i4b), dimension(neq), intent(inout)  IW,
real(dp), dimension(neq), intent(inout)  W,
real(dp), intent(in)  RELAX,
integer(i4b), intent(inout)  IPCFLAG,
real(dp), intent(in)  DELTA 
)

Update the ILU0 preconditioner using the current coefficient matrix.

Parameters
[in]njanumber of non-zero entries
[in]neqnumber of equations
[in]amatcoefficient matrix
[in]iaCRS row pointers
[in]jaCRS column pointers
[in,out]apcpreconditioned matrix
[in,out]iapcpreconditioner CRS row pointers
[in,out]japcpreconditioner CRS column pointers
[in,out]iwpreconditioner integer work vector
[in,out]wpreconditioner work vector
[in]relaxMILU0 preconditioner relaxation factor
[in,out]ipcflagpreconditioner error flag
[in]deltafactor used to correct non-diagonally dominant matrices

Definition at line 928 of file ImsLinearBase.f90.

931  ! -- dummy variables
932  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
933  integer(I4B), INTENT(IN) :: NEQ !< number of equations
934  real(DP), DIMENSION(NJA), INTENT(IN) :: AMAT !< coefficient matrix
935  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< CRS row pointers
936  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< CRS column pointers
937  real(DP), DIMENSION(NJA), INTENT(INOUT) :: APC !< preconditioned matrix
938  integer(I4B), DIMENSION(NEQ + 1), INTENT(INOUT) :: IAPC !< preconditioner CRS row pointers
939  integer(I4B), DIMENSION(NJA), INTENT(INOUT) :: JAPC !< preconditioner CRS column pointers
940  integer(I4B), DIMENSION(NEQ), INTENT(INOUT) :: IW !< preconditioner integer work vector
941  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: W !< preconditioner work vector
942  real(DP), INTENT(IN) :: RELAX !< MILU0 preconditioner relaxation factor
943  integer(I4B), INTENT(INOUT) :: IPCFLAG !< preconditioner error flag
944  real(DP), INTENT(IN) :: DELTA !< factor used to correct non-diagonally dominant matrices
945  ! -- local variables
946  integer(I4B) :: ic0, ic1
947  integer(I4B) :: iic0, iic1
948  integer(I4B) :: iu, iiu
949  integer(I4B) :: j, n
950  integer(I4B) :: jj
951  integer(I4B) :: jcol, jw
952  integer(I4B) :: jjcol
953  real(DP) :: drelax
954  real(DP) :: sd1
955  real(DP) :: tl
956  real(DP) :: rs
957  real(DP) :: d
958  !
959  ! -- initialize local variables
960  drelax = relax
961  DO n = 1, neq
962  iw(n) = 0
963  w(n) = dzero
964  END DO
965  main: DO n = 1, neq
966  ic0 = ia(n)
967  ic1 = ia(n + 1) - 1
968  DO j = ic0, ic1
969  jcol = ja(j)
970  iw(jcol) = 1
971  w(jcol) = w(jcol) + amat(j)
972  END DO
973  ic0 = iapc(n)
974  ic1 = iapc(n + 1) - 1
975  iu = japc(n)
976  rs = dzero
977  lower: DO j = ic0, iu - 1
978  jcol = japc(j)
979  iic0 = iapc(jcol)
980  iic1 = iapc(jcol + 1) - 1
981  iiu = japc(jcol)
982  tl = w(jcol) * apc(jcol)
983  w(jcol) = tl
984  DO jj = iiu, iic1
985  jjcol = japc(jj)
986  jw = iw(jjcol)
987  IF (jw .NE. 0) THEN
988  w(jjcol) = w(jjcol) - tl * apc(jj)
989  ELSE
990  rs = rs + tl * apc(jj)
991  END IF
992  END DO
993  END DO lower
994  !
995  ! -- DIAGONAL - CALCULATE INVERSE OF DIAGONAL FOR SOLUTION
996  d = w(n)
997  tl = (done + delta) * d - (drelax * rs)
998  !
999  ! -- ENSURE THAT THE SIGN OF THE DIAGONAL HAS NOT CHANGED AND IS
1000  sd1 = sign(d, tl)
1001  IF (sd1 .NE. d) THEN
1002  !
1003  ! -- USE SMALL VALUE IF DIAGONAL SCALING IS NOT EFFECTIVE FOR
1004  ! PIVOTS THAT CHANGE THE SIGN OF THE DIAGONAL
1005  IF (ipcflag > 1) THEN
1006  tl = sign(dem6, d)
1007  !
1008  ! -- DIAGONAL SCALING CONTINUES TO BE EFFECTIVE
1009  ELSE
1010  ipcflag = 1
1011  EXIT main
1012  END IF
1013  END IF
1014  IF (abs(tl) == dzero) THEN
1015  !
1016  ! -- USE SMALL VALUE IF DIAGONAL SCALING IS NOT EFFECTIVE FOR
1017  ! ZERO PIVOTS
1018  IF (ipcflag > 1) THEN
1019  tl = sign(dem6, d)
1020  !
1021  ! -- DIAGONAL SCALING CONTINUES TO BE EFFECTIVE FOR ELIMINATING
1022  ELSE
1023  ipcflag = 1
1024  EXIT main
1025  END IF
1026  END IF
1027  apc(n) = done / tl
1028  !
1029  ! -- RESET POINTER FOR IW TO ZERO
1030  iw(n) = 0
1031  w(n) = dzero
1032  DO j = ic0, ic1
1033  jcol = japc(j)
1034  apc(j) = w(jcol)
1035  iw(jcol) = 0
1036  w(jcol) = dzero
1037  END DO
1038  END DO main
1039  !
1040  ! -- RESET IPCFLAG IF SUCCESSFUL COMPLETION OF MAIN
1041  ipcflag = 0
Here is the caller graph for this function:

◆ ims_base_pcjac()

subroutine imslinearbasemodule::ims_base_pcjac ( integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NEQ,
real(dp), dimension(nja), intent(in)  AMAT,
real(dp), dimension(neq), intent(inout)  APC,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA 
)

Calculate the Jacobi preconditioner (inverse of the diagonal) using the current coefficient matrix.

Parameters
[in]njanumber of non-zero entries
[in]neqnumber of equations
[in]amatcoefficient matrix
[in,out]apcpreconditioner matrix
[in]iaCRS row pointers
[in]jaCRS column pointers

Definition at line 872 of file ImsLinearBase.f90.

873  ! -- dummy variables
874  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
875  integer(I4B), INTENT(IN) :: NEQ !< number of equations
876  real(DP), DIMENSION(NJA), INTENT(IN) :: AMAT !< coefficient matrix
877  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: APC !< preconditioner matrix
878  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< CRS row pointers
879  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< CRS column pointers
880  ! -- local variables
881  integer(I4B) :: i, n
882  integer(I4B) :: ic0, ic1
883  integer(I4B) :: id
884  real(DP) :: tv
885  ! -- code
886  DO n = 1, neq
887  ic0 = ia(n)
888  ic1 = ia(n + 1) - 1
889  id = ia(n)
890  DO i = ic0, ic1
891  IF (ja(i) == n) THEN
892  id = i
893  EXIT
894  END IF
895  END DO
896  tv = amat(id)
897  IF (abs(tv) > dzero) tv = done / tv
898  apc(n) = tv
899  END DO

◆ ims_base_pcu()

subroutine imslinearbasemodule::ims_base_pcu ( integer(i4b), intent(in)  IOUT,
integer(i4b), intent(in)  NJA,
integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NIAPC,
integer(i4b), intent(in)  NJAPC,
integer(i4b), intent(in)  IPC,
real(dp), intent(in)  RELAX,
real(dp), dimension(nja), intent(in)  AMAT,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA,
real(dp), dimension(njapc), intent(inout)  APC,
integer(i4b), dimension(niapc + 1), intent(inout)  IAPC,
integer(i4b), dimension(njapc), intent(inout)  JAPC,
integer(i4b), dimension(niapc), intent(inout)  IW,
real(dp), dimension(niapc), intent(inout)  W,
integer(i4b), intent(in)  LEVEL,
real(dp), intent(in)  DROPTOL,
integer(i4b), intent(in)  NJLU,
integer(i4b), intent(in)  NJW,
integer(i4b), intent(in)  NWLU,
integer(i4b), dimension(njlu), intent(inout)  JLU,
integer(i4b), dimension(njw), intent(inout)  JW,
real(dp), dimension(nwlu), intent(inout)  WLU 
)

Update the preconditioner using the current coefficient matrix.

Parameters
[in]ioutsimulation listing file unit
[in]njanumber of non-zero entries
[in]neqnumber of equations
[in]niapcpreconditioner number of rows
[in]njapcpreconditioner number of non-zero entries
[in]ipcprecoditioner (1) ILU0 (2) MILU0 (3) ILUT (4) MILUT
[in]relaxpreconditioner relaxation factor for MILU0 and MILUT
[in]amatcoefficient matrix
[in]iaCRS row pointers
[in]jaCRS column pointers
[in,out]apcpreconditioner matrix
[in,out]iapcpreconditioner CRS row pointers
[in,out]japcpreconditioner CRS column pointers
[in,out]iwpreconditioner integed work vector
[in,out]wpreconditioner work vector
[in]levelnumber of levels of fill for ILUT and MILUT
[in]droptoldrop tolerance
[in]njlulength of JLU working vector
[in]njwlength of JW working vector
[in]nwlulength of WLU working vector
[in,out]jluILUT/MILUT JLU working vector
[in,out]jwILUT/MILUT JW working vector
[in,out]wluILUT/MILUT WLU working vector

Definition at line 761 of file ImsLinearBase.f90.

764  ! -- modules
766  ! -- dummy variables
767  integer(I4B), INTENT(IN) :: IOUT !< simulation listing file unit
768  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
769  integer(I4B), INTENT(IN) :: NEQ !< number of equations
770  integer(I4B), INTENT(IN) :: NIAPC !< preconditioner number of rows
771  integer(I4B), INTENT(IN) :: NJAPC !< preconditioner number of non-zero entries
772  integer(I4B), INTENT(IN) :: IPC !< precoditioner (1) ILU0 (2) MILU0 (3) ILUT (4) MILUT
773  real(DP), INTENT(IN) :: RELAX !< preconditioner relaxation factor for MILU0 and MILUT
774  real(DP), DIMENSION(NJA), INTENT(IN) :: AMAT !< coefficient matrix
775  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< CRS row pointers
776  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< CRS column pointers
777  real(DP), DIMENSION(NJAPC), INTENT(INOUT) :: APC !< preconditioner matrix
778  integer(I4B), DIMENSION(NIAPC + 1), INTENT(INOUT) :: IAPC !< preconditioner CRS row pointers
779  integer(I4B), DIMENSION(NJAPC), INTENT(INOUT) :: JAPC !< preconditioner CRS column pointers
780  integer(I4B), DIMENSION(NIAPC), INTENT(INOUT) :: IW !< preconditioner integed work vector
781  real(DP), DIMENSION(NIAPC), INTENT(INOUT) :: W !< preconditioner work vector
782  ! -- ILUT dummy variables
783  integer(I4B), INTENT(IN) :: LEVEL !< number of levels of fill for ILUT and MILUT
784  real(DP), INTENT(IN) :: DROPTOL !< drop tolerance
785  integer(I4B), INTENT(IN) :: NJLU !< length of JLU working vector
786  integer(I4B), INTENT(IN) :: NJW !< length of JW working vector
787  integer(I4B), INTENT(IN) :: NWLU !< length of WLU working vector
788  integer(I4B), DIMENSION(NJLU), INTENT(INOUT) :: JLU !< ILUT/MILUT JLU working vector
789  integer(I4B), DIMENSION(NJW), INTENT(INOUT) :: JW !< ILUT/MILUT JW working vector
790  real(DP), DIMENSION(NWLU), INTENT(INOUT) :: WLU !< ILUT/MILUT WLU working vector
791  ! -- local variables
792  character(len=LINELENGTH) :: errmsg
793  character(len=100), dimension(5), parameter :: cerr = &
794  ["Elimination process has generated a row in L or U whose length is > n.", &
795  &"The matrix L overflows the array al. ", &
796  &"The matrix U overflows the array alu. ", &
797  &"Illegal value for lfil. ", &
798  &"Zero row encountered. "]
799  integer(i4b) :: ipcflag
800  integer(I4B) :: icount
801  integer(I4B) :: ierr
802  real(DP) :: delta
803  ! -- formats
804 2000 FORMAT(/, ' MATRIX IS SEVERELY NON-DIAGONALLY DOMINANT.', &
805  /, ' ADDED SMALL VALUE TO PIVOT ', i0, ' TIMES IN', &
806  ' IMSLINEARSUB_PCU.')
807  !
808  ! -- initialize local variables
809  ipcflag = 0
810  icount = 0
811  delta = dzero
812  pcscale: DO
813  SELECT CASE (ipc)
814  !
815  ! -- ILU0 AND MILU0
816  CASE (1, 2)
817  CALL ims_base_pcilu0(nja, neq, amat, ia, ja, &
818  apc, iapc, japc, iw, w, &
819  relax, ipcflag, delta)
820  !
821  ! -- ILUT AND MILUT
822  CASE (3, 4)
823  ierr = 0
824  CALL ilut(neq, amat, ja, ia, level, droptol, &
825  apc, jlu, iw, njapc, wlu, jw, ierr, &
826  relax, ipcflag, delta)
827  if (ierr /= 0) then
828  if (ierr > 0) then
829  write (errmsg, '(a,1x,i0,1x,a)') &
830  'ILUT: zero pivot encountered at step number', ierr, '.'
831  else
832  write (errmsg, '(a,1x,a)') 'ILUT:', cerr(-ierr)
833  end if
834  call store_error(errmsg)
835  call parser%StoreErrorUnit()
836  end if
837  !
838  ! -- ADDITIONAL PRECONDITIONERS
839  CASE DEFAULT
840  ipcflag = 0
841  END SELECT
842  IF (ipcflag < 1) THEN
843  EXIT pcscale
844  END IF
845  delta = 1.5d0 * delta + dem3
846  ipcflag = 0
847  IF (delta > dhalf) THEN
848  delta = dhalf
849  ipcflag = 2
850  END IF
851  icount = icount + 1
852  !
853  ! -- terminate pcscale loop if not making progress
854  if (icount > 10) then
855  exit pcscale
856  end if
857 
858  END DO pcscale
859  !
860  ! -- write error message if small value added to pivot
861  if (icount > 0) then
862  write (iout, 2000) icount
863  end if
subroutine ilut(n, a, ja, ia, lfil, droptol, alu, jlu, ju, iwk, w, jw, ierr, relax, izero, delta)
Definition: ilut.f90:51
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_residual()

subroutine imslinearbasemodule::ims_base_residual ( integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
real(dp), dimension(neq), intent(in)  X,
real(dp), dimension(neq), intent(in)  B,
real(dp), dimension(neq), intent(inout)  D,
real(dp), dimension(nja), intent(in)  A,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA 
)

Subroutine to calculate the residual.

Parameters
[in]neqlength of vectors
[in]njalength of coefficient matrix
[in]xdependent variable
[in]bright-hand side
[in,out]dresidual
[in]acoefficient matrix
[in]iaCRS row pointers
[in]jaCRS column pointers

Definition at line 1291 of file ImsLinearBase.f90.

1292  ! -- dummy variables
1293  integer(I4B), INTENT(IN) :: NEQ !< length of vectors
1294  integer(I4B), INTENT(IN) :: NJA !< length of coefficient matrix
1295  real(DP), DIMENSION(NEQ), INTENT(IN) :: X !< dependent variable
1296  real(DP), DIMENSION(NEQ), INTENT(IN) :: B !< right-hand side
1297  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: D !< residual
1298  real(DP), DIMENSION(NJA), INTENT(IN) :: A !< coefficient matrix
1299  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< CRS row pointers
1300  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< CRS column pointers
1301  ! -- local variables
1302  integer(I4B) :: n
1303  ! -- code
1304  !
1305  ! -- calculate matrix-vector product
1306  call amux(neq, x, d, a, ja, ia)
1307  !
1308  ! -- subtract matrix-vector product from right-hand side
1309  DO n = 1, neq
1310  d(n) = b(n) - d(n)
1311  END DO
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ims_base_scale()

subroutine imslinearbasemodule::ims_base_scale ( integer(i4b), intent(in)  IOPT,
integer(i4b), intent(in)  ISCL,
integer(i4b), intent(in)  NEQ,
integer(i4b), intent(in)  NJA,
integer(i4b), dimension(neq + 1), intent(in)  IA,
integer(i4b), dimension(nja), intent(in)  JA,
real(dp), dimension(nja), intent(inout)  AMAT,
real(dp), dimension(neq), intent(inout)  X,
real(dp), dimension(neq), intent(inout)  B,
real(dp), dimension(neq), intent(inout)  DSCALE,
real(dp), dimension(neq), intent(inout)  DSCALE2 
)

Scale the coefficient matrix (AMAT), the right-hand side (B), and the estimate of the dependent variable (X).

Parameters
[in]ioptflag to scale (0) or unscale the system of equations
[in]isclscaling option (1) symmetric (2) L-2 norm
[in]neqnumber of equations
[in]njanumber of non-zero entries
[in]iaCRS row pointer
[in]jaCRS column pointer
[in,out]amatcoefficient matrix
[in,out]xdependent variable
[in,out]bright-hand side
[in,out]dscalefirst scaling vector
[in,out]dscale2second scaling vector

Definition at line 619 of file ImsLinearBase.f90.

621  ! -- dummy variables
622  integer(I4B), INTENT(IN) :: IOPT !< flag to scale (0) or unscale the system of equations
623  integer(I4B), INTENT(IN) :: ISCL !< scaling option (1) symmetric (2) L-2 norm
624  integer(I4B), INTENT(IN) :: NEQ !< number of equations
625  integer(I4B), INTENT(IN) :: NJA !< number of non-zero entries
626  integer(I4B), DIMENSION(NEQ + 1), INTENT(IN) :: IA !< CRS row pointer
627  integer(I4B), DIMENSION(NJA), INTENT(IN) :: JA !< CRS column pointer
628  real(DP), DIMENSION(NJA), INTENT(INOUT) :: AMAT !< coefficient matrix
629  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: X !< dependent variable
630  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: B !< right-hand side
631  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: DSCALE !< first scaling vector
632  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: DSCALE2 !< second scaling vector
633  ! -- local variables
634  integer(I4B) :: i, n
635  integer(I4B) :: id, jc
636  integer(I4B) :: i0, i1
637  real(DP) :: v, c1, c2
638  !
639  ! -- SCALE SCALE AMAT, X, AND B
640  IF (iopt == 0) THEN
641  !
642  ! -- SYMMETRIC SCALING
643  SELECT CASE (iscl)
644  CASE (1)
645  DO n = 1, neq
646  id = ia(n)
647  v = amat(id)
648  c1 = done / sqrt(abs(v))
649  dscale(n) = c1
650  dscale2(n) = c1
651  END DO
652  !
653  ! -- SCALE AMAT -- AMAT = DSCALE(row) * AMAT(i) * DSCALE2(col)
654  DO n = 1, neq
655  c1 = dscale(n)
656  i0 = ia(n)
657  i1 = ia(n + 1) - 1
658  DO i = i0, i1
659  jc = ja(i)
660  c2 = dscale2(jc)
661  amat(i) = c1 * amat(i) * c2
662  END DO
663  END DO
664  !
665  ! -- L-2 NORM SCALING
666  CASE (2)
667  !
668  ! -- SCALE EACH ROW SO THAT THE L-2 NORM IS 1
669  DO n = 1, neq
670  c1 = dzero
671  i0 = ia(n)
672  i1 = ia(n + 1) - 1
673  DO i = i0, i1
674  c1 = c1 + amat(i) * amat(i)
675  END DO
676  c1 = sqrt(c1)
677  IF (c1 == dzero) THEN
678  c1 = done
679  ELSE
680  c1 = done / c1
681  END IF
682  dscale(n) = c1
683  !
684  ! -- INITIAL SCALING OF AMAT -- AMAT = DSCALE(row) * AMAT(i)
685  DO i = i0, i1
686  amat(i) = c1 * amat(i)
687  END DO
688  END DO
689  !
690  ! -- SCALE EACH COLUMN SO THAT THE L-2 NORM IS 1
691  DO n = 1, neq
692  dscale2(n) = dzero
693  END DO
694  c2 = dzero
695  DO n = 1, neq
696  i0 = ia(n)
697  i1 = ia(n + 1) - 1
698  DO i = i0, i1
699  jc = ja(i)
700  c2 = amat(i)
701  dscale2(jc) = dscale2(jc) + c2 * c2
702  END DO
703  END DO
704  DO n = 1, neq
705  c2 = dscale2(n)
706  IF (c2 == dzero) THEN
707  c2 = done
708  ELSE
709  c2 = done / sqrt(c2)
710  END IF
711  dscale2(n) = c2
712  END DO
713  !
714  ! -- FINAL SCALING OF AMAT -- AMAT = DSCALE2(col) * AMAT(i)
715  DO n = 1, neq
716  i0 = ia(n)
717  i1 = ia(n + 1) - 1
718  DO i = i0, i1
719  jc = ja(i)
720  c2 = dscale2(jc)
721  amat(i) = c2 * amat(i)
722  END DO
723  END DO
724  END SELECT
725  !
726  ! -- SCALE X AND B
727  DO n = 1, neq
728  c1 = dscale(n)
729  c2 = dscale2(n)
730  x(n) = x(n) / c2
731  b(n) = b(n) * c1
732  END DO
733  !
734  ! -- UNSCALE SCALE AMAT, X, AND B
735  ELSE
736  DO n = 1, neq
737  c1 = dscale(n)
738  i0 = ia(n)
739  i1 = ia(n + 1) - 1
740  !
741  ! -- UNSCALE AMAT
742  DO i = i0, i1
743  jc = ja(i)
744  c2 = dscale2(jc)
745  amat(i) = (done / c1) * amat(i) * (done / c2)
746  END DO
747  !
748  ! -- UNSCALE X AND B
749  c2 = dscale2(n)
750  x(n) = x(n) * c2
751  b(n) = b(n) / c1
752  END DO
753  END IF
Here is the caller graph for this function:

◆ ims_base_testcnvg()

subroutine imslinearbasemodule::ims_base_testcnvg ( integer(i4b), intent(in)  Icnvgopt,
integer(i4b), intent(inout)  Icnvg,
integer(i4b), intent(in)  Iiter,
real(dp), intent(in)  Dvmax,
real(dp), intent(in)  Rmax,
real(dp), intent(in)  Rmax0,
real(dp), intent(in)  Epfact,
real(dp), intent(in)  Dvclose,
real(dp), intent(in)  Rclose 
)

General routine for testing for solver convergence based on the user-specified convergence option (Icnvgopt).

Parameters
[in]icnvgoptconvergence option - see documentation for option
[in,out]icnvgflag indicating if convergence achieved (1) or not (0)
[in]iiterinner iteration number (used for strict convergence option)
[in]dvmaxmaximum dependent-variable change
[in]rmaxmaximum flow change
[in]rmax0initial flow change (initial L2-norm)
[in]epfactfactor for reducing convergence criteria in subsequent Picard iterations
[in]dvcloseMaximum depenendent-variable change allowed
[in]rcloseMaximum flow change allowed

Definition at line 1101 of file ImsLinearBase.f90.

1104  ! -- dummy variables
1105  integer(I4B), INTENT(IN) :: Icnvgopt !< convergence option - see documentation for option
1106  integer(I4B), INTENT(INOUT) :: Icnvg !< flag indicating if convergence achieved (1) or not (0)
1107  integer(I4B), INTENT(IN) :: Iiter !< inner iteration number (used for strict convergence option)
1108  real(DP), INTENT(IN) :: Dvmax !< maximum dependent-variable change
1109  real(DP), INTENT(IN) :: Rmax !< maximum flow change
1110  real(DP), INTENT(IN) :: Rmax0 !< initial flow change (initial L2-norm)
1111  real(DP), INTENT(IN) :: Epfact !< factor for reducing convergence criteria in subsequent Picard iterations
1112  real(DP), INTENT(IN) :: Dvclose !< Maximum depenendent-variable change allowed
1113  real(DP), INTENT(IN) :: Rclose !< Maximum flow change allowed
1114  ! -- code
1115  IF (icnvgopt == 0) THEN
1116  IF (abs(dvmax) <= dvclose .AND. abs(rmax) <= rclose) THEN
1117  icnvg = 1
1118  END IF
1119  ELSE IF (icnvgopt == 1) THEN
1120  IF (abs(dvmax) <= dvclose .AND. abs(rmax) <= rclose) THEN
1121  IF (iiter == 1) THEN
1122  icnvg = 1
1123  ELSE
1124  icnvg = -1
1125  END IF
1126  END IF
1127  ELSE IF (icnvgopt == 2) THEN
1128  IF (abs(dvmax) <= dvclose .OR. rmax <= rclose) THEN
1129  icnvg = 1
1130  ELSE IF (rmax <= rmax0 * epfact) THEN
1131  icnvg = -1
1132  END IF
1133  ELSE IF (icnvgopt == 3) THEN
1134  IF (abs(dvmax) <= dvclose) THEN
1135  icnvg = 1
1136  ELSE IF (rmax <= rmax0 * rclose) THEN
1137  icnvg = -1
1138  END IF
1139  ELSE IF (icnvgopt == 4) THEN
1140  IF (abs(dvmax) <= dvclose .AND. rmax <= rclose) THEN
1141  icnvg = 1
1142  ELSE IF (rmax <= rmax0 * epfact) THEN
1143  icnvg = -1
1144  END IF
1145  END IF
Here is the caller graph for this function:

◆ ims_calc_pcdims()

subroutine imslinearbasemodule::ims_calc_pcdims ( integer(i4b), intent(in)  neq,
integer(i4b), intent(in)  nja,
integer(i4b), dimension(:), intent(in)  ia,
integer(i4b), intent(in)  level,
integer(i4b), intent(in)  ipc,
integer(i4b), intent(inout)  niapc,
integer(i4b), intent(inout)  njapc,
integer(i4b), intent(inout)  njlu,
integer(i4b), intent(inout)  njw,
integer(i4b), intent(inout)  nwlu 
)
Parameters
[in]neqnr. of rows A
[in]njanr. of nonzeros A
[in]iaCSR row pointers A
[in]levelfill level ILU
[in]ipcIMS preconditioner type
[in,out]niapcwork array size
[in,out]njapcwork array size
[in,out]njluwork array size
[in,out]njwwork array size
[in,out]nwluwork array size

Definition at line 1148 of file ImsLinearBase.f90.

1150  integer(I4B), intent(in) :: neq !< nr. of rows A
1151  integer(I4B), intent(in) :: nja !< nr. of nonzeros A
1152  integer(I4B), dimension(:), intent(in) :: ia !< CSR row pointers A
1153  integer(I4B), intent(in) :: level !< fill level ILU
1154  integer(I4B), intent(in) :: ipc !< IMS preconditioner type
1155  integer(I4B), intent(inout) :: niapc !< work array size
1156  integer(I4B), intent(inout) :: njapc !< work array size
1157  integer(I4B), intent(inout) :: njlu !< work array size
1158  integer(I4B), intent(inout) :: njw !< work array size
1159  integer(I4B), intent(inout) :: nwlu !< work array size
1160  ! local
1161  integer(I4B) :: n, i
1162  integer(I4B) :: ijlu, ijw, iwlu, iwk
1163 
1164  ijlu = 1
1165  ijw = 1
1166  iwlu = 1
1167 
1168  ! ILU0 and MILU0
1169  niapc = neq
1170  njapc = nja
1171 
1172  ! ILUT and MILUT
1173  if (ipc == 3 .or. ipc == 4) then
1174  niapc = neq
1175  if (level > 0) then
1176  iwk = neq * (level * 2 + 1)
1177  else
1178  iwk = 0
1179  do n = 1, neq
1180  i = ia(n + 1) - ia(n)
1181  if (i > iwk) then
1182  iwk = i
1183  end if
1184  end do
1185  iwk = neq * iwk
1186  end if
1187  njapc = iwk
1188  ijlu = iwk
1189  ijw = 2 * neq
1190  iwlu = neq + 1
1191  end if
1192 
1193  njlu = ijlu
1194  njw = ijw
1195  nwlu = iwlu
1196 
Here is the caller graph for this function:

Variable Documentation

◆ parser

type(blockparsertype), private imslinearbasemodule::parser
private

Definition at line 19 of file ImsLinearBase.f90.

19  type(BlockParserType), private :: parser