MODFLOW 6  version 6.6.0.dev0
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...
 

Function/Subroutine Documentation

◆ 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 17 of file ImsLinearMisc.f90.

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