-
Notifications
You must be signed in to change notification settings - Fork 2
/
Copy pathtridiagonal_cyclic.f90
276 lines (276 loc) · 11.5 KB
/
tridiagonal_cyclic.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
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
!
! SCID-TDSE: Simple 1-electron atomic TDSE solver
! Copyright (C) 2015-2021 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/>.
!
!
! Solution of tridiagonal linear systems using unpivoted cyclic-reduction algorthm.
!
! See Gander and Golub, "Cyclic Reduction - History and Applications", Proceedings of
! the Workshop on Scientific Computing 10-12 March 1997 (Springer, NY), available at
! http://www.inf.ethz.ch/personal/gander/papers/cyclic.pdf
!
! We choose bot to provide the *_x (high-precision on low-precision build) version
! of this solver. Numerical properties of the cyclic-reduction solver are the
! same as for the unpivoted solver, so there is no point in implementing it for
! calls not on the critical path.
!
module tridiagonal_cyclic
use accuracy
use constants
use timer
implicit none
private
public m3c_decompose, m3c_solve
public rcsid_tridiagonal_cyclic
!
character(len=clen), save :: rcsid_tridiagonal_cyclic = "$Id: tridiagonal_cyclic.f90,v 1.2 2021/04/26 15:44:44 ps Exp $"
!
integer(ik), parameter :: max_cr_steps = 64 ! Enough to exhaust address space on a 64-bit system several times over
!
! Public interfaces
!
interface m3c_decompose
module procedure m3c_decompose_r
module procedure m3c_decompose_c
end interface m3c_decompose
!
interface m3c_solve
module procedure m3c_solve_rr
module procedure m3c_solve_rc
module procedure m3c_solve_cc
module procedure m3c_solve_rr_g
module procedure m3c_solve_rc_g
module procedure m3c_solve_cc_g
end interface m3c_solve
!
! Internal interfaces
!
interface cr_reduce_step
module procedure cr_reduce_step_r
module procedure cr_reduce_step_c
end interface cr_reduce_step
!
interface cr_solve_forward_step
module procedure cr_solve_forward_step_rr
module procedure cr_solve_forward_step_rc
module procedure cr_solve_forward_step_cc
module procedure cr_solve_forward_step_rr_g
module procedure cr_solve_forward_step_rc_g
module procedure cr_solve_forward_step_cc_g
end interface cr_solve_forward_step
!
interface cr_solve_backward_step
module procedure cr_solve_backward_step_rr
module procedure cr_solve_backward_step_rc
module procedure cr_solve_backward_step_cc
module procedure cr_solve_backward_step_rr_g
module procedure cr_solve_backward_step_rc_g
module procedure cr_solve_backward_step_cc_g
end interface cr_solve_backward_step
!
contains
!
subroutine m3c_decompose_r(m,mf,fail)
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix, assumed to be diagonally dominant
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.
include 'tridiagonal_cyclic_m3c_decompose_common.f90'
end subroutine m3c_decompose_r
!
subroutine m3c_decompose_c(m,mf,fail)
complex(rk), intent(in) :: m (:,:) ! Tridiagonal matrix, assumed to be diagonally dominant
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.
include 'tridiagonal_cyclic_m3c_decompose_common.f90'
end subroutine m3c_decompose_c
!
subroutine m3c_solve_rr(mf,rhs,x,scr)
real(rk), intent(in) :: mf(:,:) ! Factorization data
real(rk), intent(in) :: rhs(:) ! Right-hand side
real(rk), intent(out) :: x(:) ! Solution of the linear system
real(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_common.f90'
end subroutine m3c_solve_rr
!
subroutine m3c_solve_rc(mf,rhs,x,scr)
real(rk), intent(in) :: mf(:,:) ! Factorization data
complex(rk), intent(in) :: rhs(:) ! Right-hand side
complex(rk), intent(out) :: x(:) ! Solution of the linear system
complex(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_common.f90'
end subroutine m3c_solve_rc
!
subroutine m3c_solve_cc(mf,rhs,x,scr)
complex(rk), intent(in) :: mf(:,:) ! Factorization data
complex(rk), intent(in) :: rhs(:) ! Right-hand side
complex(rk), intent(out) :: x(:) ! Solution of the linear system
complex(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_common.f90'
end subroutine m3c_solve_cc
!
subroutine m3c_solve_rr_g(mf,rhs,x,scr)
real(rk), intent(in) :: mf(:,:) ! Factorization data
real(rk), intent(in) :: rhs(:,:) ! Right-hand sides
real(rk), intent(out) :: x(:,:) ! Solutions of the linear system
real(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_g_common.f90'
end subroutine m3c_solve_rr_g
!
subroutine m3c_solve_rc_g(mf,rhs,x,scr)
real(rk), intent(in) :: mf(:,:) ! Factorization data
complex(rk), intent(in) :: rhs(:,:) ! Right-hand sides
complex(rk), intent(out) :: x(:,:) ! Solutions of the linear system
complex(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_g_common.f90'
end subroutine m3c_solve_rc_g
!
subroutine m3c_solve_cc_g(mf,rhs,x,scr)
complex(rk), intent(in) :: mf(:,:) ! Factorization data
complex(rk), intent(in) :: rhs(:,:) ! Right-hand sides
complex(rk), intent(out) :: x(:,:) ! Solutions of the linear system
complex(rk), intent(out) :: scr(:,:) ! Scratch space
!
include 'tridiagonal_cyclic_m3c_solve_g_common.f90'
end subroutine m3c_solve_cc_g
!
! Internal subroutines below this point
!
subroutine cr_reduce_step_r(m,f,fail)
real(rk), intent(in) :: m (:,:) ! Tridiagonal matrix to be reduced
real(rk), intent(out) :: f(:,:) ! Reduction data.
logical, intent(out) :: fail
!
include 'tridiagonal_cyclic_cr_reduce_step_common.f90'
end subroutine cr_reduce_step_r
!
subroutine cr_reduce_step_c(m,f,fail)
complex(rk), intent(in) :: m (:,:) ! Tridiagonal matrix to be reduced
complex(rk), intent(out) :: f(:,:) ! Reduction data.
logical, intent(out) :: fail
!
include 'tridiagonal_cyclic_cr_reduce_step_common.f90'
end subroutine cr_reduce_step_c
!
subroutine cr_solve_forward_step_rr(f,r,rr)
real(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
real(rk), intent(in) :: r (:) ! Full RHS
real(rk), intent(out) :: rr(:) ! Reduced RHS on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_common.f90'
end subroutine cr_solve_forward_step_rr
!
subroutine cr_solve_forward_step_rc(f,r,rr)
real(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:) ! Full RHS
complex(rk), intent(out) :: rr(:) ! Reduced RHS on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_common.f90'
end subroutine cr_solve_forward_step_rc
!
subroutine cr_solve_forward_step_cc(f,r,rr)
complex(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:) ! Full RHS
complex(rk), intent(out) :: rr(:) ! Reduced RHS on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_common.f90'
end subroutine cr_solve_forward_step_cc
!
subroutine cr_solve_forward_step_rr_g(f,r,rr)
real(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
real(rk), intent(in) :: r (:,:) ! Full RHS set
real(rk), intent(out) :: rr(:,:) ! Reduced RHS set on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_g_common.f90'
end subroutine cr_solve_forward_step_rr_g
!
subroutine cr_solve_forward_step_rc_g(f,r,rr)
real(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:,:) ! Full RHS set
complex(rk), intent(out) :: rr(:,:) ! Reduced RHS set on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_g_common.f90'
end subroutine cr_solve_forward_step_rc_g
!
subroutine cr_solve_forward_step_cc_g(f,r,rr)
complex(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:,:) ! Full RHS set
complex(rk), intent(out) :: rr(:,:) ! Reduced RHS set on this step
!
include 'tridiagonal_cyclic_cr_solve_forward_g_common.f90'
end subroutine cr_solve_forward_step_cc_g
!
subroutine cr_solve_backward_step_rr(f,r,xr,x)
real(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
real(rk), intent(in) :: r (:) ! RHS
real(rk), intent(in) :: xr(:) ! Reduced solution vector
real(rk), intent(out) :: x (:) ! Full Solution vector
!
include 'tridiagonal_cyclic_cr_solve_backward_common.f90'
end subroutine cr_solve_backward_step_rr
!
subroutine cr_solve_backward_step_rc(f,r,xr,x)
real(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:) ! RHS
complex(rk), intent(in) :: xr(:) ! Reduced solution vector
complex(rk), intent(out) :: x (:) ! Full Solution vector
!
include 'tridiagonal_cyclic_cr_solve_backward_common.f90'
end subroutine cr_solve_backward_step_rc
!
subroutine cr_solve_backward_step_cc(f,r,xr,x)
complex(rk), intent(in) :: f(:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:) ! RHS
complex(rk), intent(in) :: xr(:) ! Reduced solution vector
complex(rk), intent(out) :: x (:) ! Full Solution vector
!
include 'tridiagonal_cyclic_cr_solve_backward_common.f90'
end subroutine cr_solve_backward_step_cc
!
subroutine cr_solve_backward_step_rr_g(f,r,xr,x)
real(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
real(rk), intent(in) :: r (:,:) ! RHS set
real(rk), intent(in) :: xr(:,:) ! Reduced solution vectors
real(rk), intent(out) :: x (:,:) ! Full Solution vectors
!
include 'tridiagonal_cyclic_cr_solve_backward_g_common.f90'
end subroutine cr_solve_backward_step_rr_g
!
subroutine cr_solve_backward_step_rc_g(f,r,xr,x)
real(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:,:) ! RHS set
complex(rk), intent(in) :: xr(:,:) ! Reduced solution vectors
complex(rk), intent(out) :: x (:,:) ! Full Solution vectors
!
include 'tridiagonal_cyclic_cr_solve_backward_g_common.f90'
end subroutine cr_solve_backward_step_rc_g
!
subroutine cr_solve_backward_step_cc_g(f,r,xr,x)
complex(rk), intent(in) :: f (:,:) ! Reduction data from cr_reduce_step
complex(rk), intent(in) :: r (:,:) ! RHS set
complex(rk), intent(in) :: xr(:,:) ! Reduced solution vectors
complex(rk), intent(out) :: x (:,:) ! Full Solution vectors
!
include 'tridiagonal_cyclic_cr_solve_backward_g_common.f90'
end subroutine cr_solve_backward_step_cc_g
!
end module tridiagonal_cyclic