-
Notifications
You must be signed in to change notification settings - Fork 0
/
tridiagonal_tools.f90
260 lines (260 loc) · 9.77 KB
/
tridiagonal_tools.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
!
! hydrogen-tunnel: Static-field tunneling in a central potential
! Copyright (C) 2018-2022 Serguei Patchkovskii, [email protected]
!
! This program is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! This program is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with this program. If not, see <https://www.gnu.org/licenses/>.
!
! Simple tools for dealing with three-diagonal matrices and systems of equations.
!
! All tridiagonal matrices used by this module are in the following format:
! m(1,i) = Main diagonal; element (i,i)
! m(2,i) = Sub-diagonal; element (i+1,i)
! m(3,i) = Super-diagonal; element (i,i+1)
!
module tridiagonal_tools
use accuracy
use constants
use timer
implicit none
private
public m3d_multiply, m3d_decompose, m3d_solve
public m3d_multiply_x, m3d_decompose_x, m3d_solve_x
public m3d_left_scale, m3d_right_scale, m3d_transpose
public rcsid_tridiagonal_tools
!
character(len=clen), save :: rcsid_tridiagonal_tools = "$Id: $"
!
interface m3d_multiply
module procedure m3d_multiply_rr
module procedure m3d_multiply_rc
module procedure m3d_multiply_cc
module procedure m3d_multiply_rr_g
module procedure m3d_multiply_rc_g
module procedure m3d_multiply_cc_g
end interface m3d_multiply
interface m3d_multiply_x
module procedure m3d_multiply_xrr_g
end interface m3d_multiply_x
!
interface m3d_decompose
module procedure m3d_decompose_r
module procedure m3d_decompose_c
end interface m3d_decompose
interface m3d_decompose_x
module procedure m3d_decompose_xr
end interface m3d_decompose_x
!
interface m3d_solve
module procedure m3d_solve_rr
module procedure m3d_solve_rc
module procedure m3d_solve_cc
module procedure m3d_solve_rr_g
module procedure m3d_solve_rc_g
module procedure m3d_solve_cc_g
end interface m3d_solve
interface m3d_solve_x
module procedure m3d_solve_xrr_g
end interface m3d_solve_x
!
interface m3d_left_scale
module procedure m3d_left_scale_cr
end interface m3d_left_scale
!
interface m3d_right_scale
module procedure m3d_right_scale_rc
end interface m3d_right_scale
!
interface m3d_transpose
module procedure m3d_transpose_r
end interface m3d_transpose
!
contains
!
subroutine m3d_multiply_rr(m,v,mv)
real(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
real(rk), intent(in) :: v(:) ! Vector
real(rk), intent(out) :: mv(:) ! Vector
!
include "tridiagonal_tools_m3d_multiply_common.f90"
end subroutine m3d_multiply_rr
!
subroutine m3d_multiply_rc(m,v,mv)
real(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
complex(rk), intent(in) :: v(:) ! Vector
complex(rk), intent(out) :: mv(:) ! Vector
!
include "tridiagonal_tools_m3d_multiply_common.f90"
end subroutine m3d_multiply_rc
!
subroutine m3d_multiply_cc(m,v,mv)
complex(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
complex(rk), intent(in) :: v(:) ! Vector
complex(rk), intent(out) :: mv(:) ! Vector
!
include "tridiagonal_tools_m3d_multiply_common.f90"
end subroutine m3d_multiply_cc
!
subroutine m3d_multiply_rr_g(m,v,mv)
real(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
real(rk), intent(in) :: v(:,:) ! Vectors
real(rk), intent(out) :: mv(:,:) ! Vectors
!
include "tridiagonal_tools_m3d_multiply_g_common.f90"
end subroutine m3d_multiply_rr_g
!
subroutine m3d_multiply_xrr_g(m,v,mv)
real(xk), intent(in) :: m(:,:) ! Tri-diagonal matrix
real(xk), intent(in) :: v(:,:) ! Vectors
real(xk), intent(out) :: mv(:,:) ! Vectors
!
include "tridiagonal_tools_m3d_multiply_g_common.f90"
end subroutine m3d_multiply_xrr_g
!
subroutine m3d_multiply_rc_g(m,v,mv)
real(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
complex(rk), intent(in) :: v(:,:) ! Vectors
complex(rk), intent(out) :: mv(:,:) ! Vectors
!
include "tridiagonal_tools_m3d_multiply_g_common.f90"
end subroutine m3d_multiply_rc_g
!
subroutine m3d_multiply_cc_g(m,v,mv)
complex(rk), intent(in) :: m(:,:) ! Tri-diagonal matrix
complex(rk), intent(in) :: v(:,:) ! Vectors
complex(rk), intent(out) :: mv(:,:) ! Vectors
!
include "tridiagonal_tools_m3d_multiply_g_common.f90"
end subroutine m3d_multiply_cc_g
!
subroutine m3d_decompose_r(m,mf,fail)
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix, assumed to be diagonally dominant
! As the result, we do not need to bother with pivoting
real(rk), intent(out) :: mf(:,:) ! Factorized tridiagonal matrix
logical, intent(out), optional :: fail ! Set to .true. if decomposition fails;
! if fail is absent, abort on decomposition failure.
!
real(rk) :: denom
!
include "tridiagonal_tools_m3d_decompose_common.f90"
end subroutine m3d_decompose_r
!
subroutine m3d_decompose_xr(m,mf,fail)
real(xk), intent(in) :: m (:,:) ! Tridiagonal matrix, assumed to be diagonally dominant
! As the result, we do not need to bother with pivoting
real(xk), intent(out) :: mf(:,:) ! Factorized tridiagonal matrix
logical, intent(out), optional :: fail ! Set to .true. if decomposition fails;
! if fail is absent, abort on decomposition failure.
!
real(xk) :: denom
!
include "tridiagonal_tools_m3d_decompose_common.f90"
end subroutine m3d_decompose_xr
!
subroutine m3d_decompose_c(m,mf,fail)
complex(rk), intent(in) :: m (:,:) ! Tridiagonal matrix, assumed to be diagonally dominant
! As the result, we do not need to bother with pivoting
complex(rk), intent(out) :: mf(:,:) ! Factorized tridiagonal matrix
logical, intent(out), optional :: fail ! Set to .true. if decomposition fails;
! if fail is absent, abort on decomposition failure.
!
complex(rk) :: denom
!
include "tridiagonal_tools_m3d_decompose_common.f90"
end subroutine m3d_decompose_c
!
subroutine m3d_solve_rr(mf,r,x)
real(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
real(rk), intent(in) :: r (:) ! Right-hand size
real(rk), intent(out) :: x (:) ! Solution vector
!
include "tridiagonal_tools_m3d_solve_common.f90"
end subroutine m3d_solve_rr
!
subroutine m3d_solve_rc(mf,r,x)
real(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
complex(rk), intent(in) :: r (:) ! Right-hand size
complex(rk), intent(out) :: x (:) ! Solution vector
!
include "tridiagonal_tools_m3d_solve_common.f90"
end subroutine m3d_solve_rc
!
subroutine m3d_solve_cc(mf,r,x)
complex(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
complex(rk), intent(in) :: r (:) ! Right-hand size
complex(rk), intent(out) :: x (:) ! Solution vector
!
include "tridiagonal_tools_m3d_solve_common.f90"
end subroutine m3d_solve_cc
!
subroutine m3d_solve_rr_g(mf,r,x)
real(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
real(rk), intent(in) :: r (:,:) ! Right-hand sizes
real(rk), intent(out) :: x (:,:) ! Solution vectors
!
include "tridiagonal_tools_m3d_solve_g_common.f90"
end subroutine m3d_solve_rr_g
!
subroutine m3d_solve_xrr_g(mf,r,x)
real(xk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
real(xk), intent(in) :: r (:,:) ! Right-hand sizes
real(xk), intent(out) :: x (:,:) ! Solution vectors
!
include "tridiagonal_tools_m3d_solve_g_common.f90"
end subroutine m3d_solve_xrr_g
!
subroutine m3d_solve_rc_g(mf,r,x)
real(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
complex(rk), intent(in) :: r (:,:) ! Right-hand sizes
complex(rk), intent(out) :: x (:,:) ! Solution vectors
!
include "tridiagonal_tools_m3d_solve_g_common.f90"
end subroutine m3d_solve_rc_g
!
subroutine m3d_solve_cc_g(mf,r,x)
complex(rk), intent(in) :: mf(:,:) ! Factorized tridiagonal matrix
complex(rk), intent(in) :: r (:,:) ! Right-hand sizes
complex(rk), intent(out) :: x (:,:) ! Solution vectors
!
include "tridiagonal_tools_m3d_solve_g_common.f90"
end subroutine m3d_solve_cc_g
!
subroutine m3d_left_scale_cr(s,m,sm)
complex(rk), intent(in) :: s (:) ! Diagonal matrix
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix
complex(rk), intent(out) :: sm(:,:) ! Tridiagonal s . m
!
include "tridiagonal_tools_m3d_left_scale_common.f90"
end subroutine m3d_left_scale_cr
!
subroutine m3d_right_scale_rc(m,s,ms)
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix
complex(rk), intent(in) :: s (:) ! Diagonal matrix
complex(rk), intent(out) :: ms(:,:) ! Tridiagonal m . s
!
include "tridiagonal_tools_m3d_right_scale_common.f90"
end subroutine m3d_right_scale_rc
!
subroutine m3d_transpose_r(m,mt)
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix
real(rk), intent(out) :: mt(:,:) ! Transpose of m
!
if (any(ubound(m)/=ubound(mt))) then
stop 'tridiagonal_tools%m3d_transpose_r - bad input sizes'
end if
!
mt(1,:) = m(1,:)
mt(2,:) = m(3,:)
mt(3,:) = m(2,:)
end subroutine m3d_transpose_r
end module tridiagonal_tools