MODFLOW 6  version 6.7.0.dev0
USGS Modular Hydrologic Model
ImsLinearMisc.f90
Go to the documentation of this file.
2 
3  use kindmodule, only: dp, i4b
4  use constantsmodule, only: dzero, done
5 
6  private
7  public :: ims_misc_thomas
8 
9 CONTAINS
10 
11  !> @brief Tridiagonal solve using the Thomas algorithm
12  !!
13  !! Subroutine to solve tridiagonal linear equations using the
14  !! Thomas algorithm.
15  !!
16  !<
17  subroutine ims_misc_thomas(n, tl, td, tu, b, x, w)
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
50  end subroutine ims_misc_thomas
51 
52 END MODULE imslinearmisc
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:65
real(dp), parameter done
real constant 1
Definition: Constants.f90:76
subroutine, public ims_misc_thomas(n, tl, td, tu, b, x, w)
Tridiagonal solve using the Thomas algorithm.
This module defines variable data types.
Definition: kind.f90:8