SCALE-RM
scale_matrix.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
11  !-----------------------------------------------------------------------------
12  !
13  !++ used modules
14  !
15  use scale_precision
16  use scale_stdio
17  !-----------------------------------------------------------------------------
18  implicit none
19  private
20  !-----------------------------------------------------------------------------
21  !
22  !++ Public procedure
23  !
24  public :: matrix_solver_tridiagonal
25 
26  interface matrix_solver_tridiagonal
27  module procedure matrix_solver_tridiagonal_1d
28  module procedure matrix_solver_tridiagonal_3d
29  end interface matrix_solver_tridiagonal
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35  !-----------------------------------------------------------------------------
36  !
37  !++ Private procedure
38  !
39  !-----------------------------------------------------------------------------
40  !
41  !++ Private parameters & variables
42  !
43  !-----------------------------------------------------------------------------
44 contains
45  !-----------------------------------------------------------------------------
47  subroutine matrix_solver_tridiagonal_1d( &
48  KA, &
49  ud, &
50  md, &
51  ld, &
52  iv, &
53  ov )
54  implicit none
55 
56  integer, intent(in) :: ka ! array size
57  real(RP), intent(in) :: ud(ka) ! upper diagonal
58  real(RP), intent(in) :: md(ka) ! middle diagonal
59  real(RP), intent(in) :: ld(ka) ! lower diagonal
60  real(RP), intent(in) :: iv(ka) ! input vector
61  real(RP), intent(out) :: ov(ka) ! output vector
62 
63  real(RP) :: c(ka)
64  real(RP) :: d(ka)
65 
66  integer :: k
67  !---------------------------------------------------------------------------
68 
69  ! foward reduction
70  c(1) = ud(1) / md(1)
71  d(1) = iv(1) / md(1)
72  do k = 2, ka
73  c(k) = ud(k) / ( md(k) - ld(k) * c(k-1) )
74  d(k) = ( iv(k) - ld(k) * d(k-1) ) / ( md(k) - ld(k) * c(k-1) )
75  enddo
76 
77  ! backward substitution
78  ov(ka) = d(ka)
79  do k = ka-1, 1, -1
80  ov(k) = d(k) - c(k) * ov(k+1)
81  enddo
82 
83  return
84  end subroutine matrix_solver_tridiagonal_1d
85 
86  !-----------------------------------------------------------------------------
88  subroutine matrix_solver_tridiagonal_3d( &
89  KA, &
90  IA, IS, IE, &
91  JA, JS, JE, &
92  ud, &
93  md, &
94  ld, &
95  iv, &
96  ov )
97  implicit none
98 
99  integer, intent(in) :: ka ! array size
100  integer, intent(in) :: ia, is, ie ! array size
101  integer, intent(in) :: ja, js, je ! array size
102  real(RP), intent(in) :: ud(ka,ia,ja) ! upper diagonal
103  real(RP), intent(in) :: md(ka,ia,ja) ! middle diagonal
104  real(RP), intent(in) :: ld(ka,ia,ja) ! lower diagonal
105  real(RP), intent(in) :: iv(ka,ia,ja) ! input vector
106  real(RP), intent(out) :: ov(ka,ia,ja) ! output vector
107 
108  real(RP) :: c(ka,ia,ja)
109  real(RP) :: d(ka,ia,ja)
110 
111  integer :: k, i, j
112  !---------------------------------------------------------------------------
113 
114  ! foward reduction
115  do j = js, je
116  do i = is, ie
117  c(1,i,j) = ud(1,i,j) / md(1,i,j)
118  d(1,i,j) = iv(1,i,j) / md(1,i,j)
119  enddo
120  enddo
121 
122  do j = js, je
123  do i = is, ie
124  do k = 2, ka
125  c(k,i,j) = ud(k,i,j) &
126  / ( md(k,i,j) - ld(k,i,j) * c(k-1,i,j) )
127  d(k,i,j) = ( iv(k,i,j) - ld(k,i,j) * d(k-1,i,j) ) &
128  / ( md(k,i,j) - ld(k,i,j) * c(k-1,i,j) )
129  enddo
130  enddo
131  enddo
132 
133  ! backward substitution
134  do j = js, je
135  do i = is, ie
136  ov(ka,i,j) = d(ka,i,j)
137  enddo
138  enddo
139 
140  do j = js, je
141  do i = is, ie
142  do k = ka-1, 1, -1
143  ov(k,i,j) = d(k,i,j) - c(k,i,j) * ov(k+1,i,j)
144  enddo
145  enddo
146  enddo
147 
148  return
149  end subroutine matrix_solver_tridiagonal_3d
150 
151 end module scale_matrix
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
module STDIO
Definition: scale_stdio.F90:12
integer, public ia
of x whole cells (local, with HALO)
integer, public ka
of z whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
module MATRIX
integer, public ie
end point of inner domain: x, local
module PRECISION
integer, public ja
of y whole cells (local, with HALO)