MODFLOW 6  version 6.7.0.dev2
USGS Modular Hydrologic Model
imslinearmisc Module Reference

Functions/Subroutines

subroutine, public ims_misc_thomas (n, tl, td, tu, b, x, w)
 Tridiagonal solve using the Thomas algorithm. More...
 
subroutine, public ims_misc_dvscale (IOPT, NEQ, DSCALE, X, B)
 @ brief Scale X and RHS More...
 

Function/Subroutine Documentation

◆ ims_misc_dvscale()

subroutine, public imslinearmisc::ims_misc_dvscale ( integer(i4b), intent(in)  IOPT,
integer(i4b), intent(in)  NEQ,
real(dp), intent(inout)  DSCALE,
real(dp), dimension(neq), intent(inout)  X,
real(dp), dimension(neq), intent(inout)  B 
)

Scale X and B to avoid big or small values. Scaling value is the maximum ABS(X).

Parameters
[in]ioptflag to scale (0) or unscale the system of equations
[in]neqnumber of equations
[in,out]dscalescaling value
[in,out]xdependent variable
[in,out]bright-hand side

Definition at line 59 of file ImsLinearMisc.f90.

60  ! -- dummy variables
61  integer(I4B), INTENT(IN) :: IOPT !< flag to scale (0) or unscale the system of equations
62  integer(I4B), INTENT(IN) :: NEQ !< number of equations
63  real(DP), INTENT(INOUT) :: DSCALE !< scaling value
64  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: X !< dependent variable
65  real(DP), DIMENSION(NEQ), INTENT(INOUT) :: B !< right-hand side
66  ! -- local variables
67  integer(I4B) :: n
68  !
69  ! -- SCALE SCALE AMAT, X, AND B
70  IF (iopt == 0) THEN
71  dscale = abs(dscale)
72  if (dscale < dem30) then
73  dscale = done
74  end if
75  DO n = 1, neq
76  b(n) = b(n) / dscale
77  x(n) = x(n) / dscale
78  END DO
79  !
80  ! -- UNSCALE X, AND B
81  ELSE
82  DO n = 1, neq
83  b(n) = b(n) * dscale
84  x(n) = x(n) * dscale
85  END DO
86  END IF
Here is the caller graph for this function:

◆ ims_misc_thomas()

subroutine, public imslinearmisc::ims_misc_thomas ( integer(i4b), intent(in)  n,
real(dp), dimension(n), intent(in)  tl,
real(dp), dimension(n), intent(in)  td,
real(dp), dimension(n), intent(in)  tu,
real(dp), dimension(n), intent(in)  b,
real(dp), dimension(n), intent(inout)  x,
real(dp), dimension(n), intent(inout)  w 
)

Subroutine to solve tridiagonal linear equations using the Thomas algorithm.

Parameters
[in]nnumber of matrix rows
[in]tllower matrix terms
[in]tddiagonal matrix terms
[in]tuupper matrix terms
[in]bright-hand side vector
[in,out]xsolution vector
[in,out]wwork vector

Definition at line 18 of file ImsLinearMisc.f90.

19  implicit none
20  ! -- dummy variables
21  integer(I4B), intent(in) :: n !< number of matrix rows
22  real(DP), dimension(n), intent(in) :: tl !< lower matrix terms
23  real(DP), dimension(n), intent(in) :: td !< diagonal matrix terms
24  real(DP), dimension(n), intent(in) :: tu !< upper matrix terms
25  real(DP), dimension(n), intent(in) :: b !< right-hand side vector
26  real(DP), dimension(n), intent(inout) :: x !< solution vector
27  real(DP), dimension(n), intent(inout) :: w !< work vector
28  ! -- local variables
29  integer(I4B) :: j
30  real(DP) :: bet
31  real(DP) :: beti
32  !
33  ! -- initialize variables
34  w(1) = dzero
35  bet = td(1)
36  beti = done / bet
37  x(1) = b(1) * beti
38  !
39  ! -- decomposition and forward substitution
40  do j = 2, n
41  w(j) = tu(j - 1) * beti
42  bet = td(j) - tl(j) * w(j)
43  beti = done / bet
44  x(j) = (b(j) - tl(j) * x(j - 1)) * beti
45  end do
46  !
47  ! -- backsubstitution
48  do j = n - 1, 1, -1
49  x(j) = x(j) - w(j + 1) * x(j + 1)
50  end do
Here is the caller graph for this function: