49 real(kind=
rp),
allocatable :: d(:)
50 real(kind=
rp),
allocatable :: w(:)
51 real(kind=
rp),
allocatable :: r(:)
52 real(kind=
rp) :: omega
55 integer :: max_iter = 10
56 logical :: recompute_diag = .true.
65 real(kind=
rp),
allocatable :: d(:)
66 real(kind=
rp),
allocatable :: w(:)
67 real(kind=
rp),
allocatable :: r(:)
68 real(kind=
rp) :: tha, dlt
71 integer :: power_its = 250
72 integer :: max_iter = 10
73 logical :: recompute_eigs = .true.
88 integer,
intent(in) :: n
89 integer,
intent(in) :: lvl
90 integer,
intent(in) :: max_iter
97 this%max_iter = max_iter
98 this%recompute_eigs = .true.
111 integer,
intent(in) :: n
112 real(kind=
rp) :: lam, b, a, rn
113 real(kind=
rp),
parameter :: boost = 1.1_rp
114 real(kind=
rp),
parameter :: lam_factor = 30.0_rp
115 real(kind=
rp) :: wtw, dtw, dtd
117 associate(w => this%w, d => this%d, coef => amg%coef, gs_h => amg%gs_h, &
118 msh=>amg%msh, xh=>amg%Xh, blst=>amg%blst)
126 if (this%lvl .eq. 0)
then
127 call gs_h%op(d, n, gs_op_add)
128 call blst%apply(d, n)
132 do i = 1, this%power_its
134 call amg%matvec(w, d, this%lvl)
136 if (this%lvl .eq. 0)
then
138 wtw =
glsc3(w, coef%mult, w, n)
143 call cmult2(d, w, 1.0_rp/sqrt(wtw), n)
145 if (this%lvl .eq. 0)
then
151 call amg%matvec(w, d, this%lvl)
152 if (this%lvl .eq. 0)
then
156 if (this%lvl .eq. 0)
then
157 dtw =
glsc3(d, coef%mult, w, n)
158 dtd =
glsc3(d, coef%mult, d, n)
166 this%tha = (b+a)/2.0_rp
167 this%dlt = (b-a)/2.0_rp
169 this%recompute_eigs = .false.
183 integer,
intent(in) :: n
184 real(kind=rp),
dimension(n),
intent(inout) :: x
185 real(kind=rp),
dimension(n),
intent(inout) :: f
186 class(tamg_hierarchy_t),
intent(inout) :: amg
187 type(ksp_monitor_t) :: ksp_results
188 integer,
optional,
intent(in) :: niter
189 integer :: iter, max_iter
190 real(kind=rp) :: rtr, rnorm
191 real(kind=rp) :: rhok, rhokp1, s1, thet, delt, tmp1, tmp2
193 if (this%recompute_eigs)
then
194 call this%comp_eig(amg, n)
197 if (
present(niter))
then
200 max_iter = this%max_iter
203 associate( w => this%w, r => this%r, d => this%d, blst=>amg%blst)
206 call amg%matvec(w, x, this%lvl)
215 call cmult(d, 1.0_rp/thet, n)
217 do iter = 1, max_iter
221 call amg%matvec(w, d, this%lvl)
224 rhokp1 = 1.0_rp / (2.0_rp * s1 - rhok)
226 tmp2 = 2.0_rp * rhokp1 / delt
229 call cmult(d, tmp1, n)
230 call add2s2(d, r, tmp2, n)
242 integer,
intent(in) :: n
243 integer,
intent(in) :: lvl
244 integer,
intent(in) :: max_iter
251 this%max_iter = max_iter
261 type(tamg_hierarchy_t),
intent(inout) :: amg
262 integer,
intent(in) :: n
266 call tamg_sample_matrix_val(val, amg, this%lvl, i, i)
267 this%d(i) = 1.0_rp / val
269 this%recompute_diag = .false.
279 integer,
intent(in) :: n
280 real(kind=rp),
dimension(n),
intent(inout) :: x
281 real(kind=rp),
dimension(n),
intent(inout) :: f
282 class(tamg_hierarchy_t),
intent(inout) :: amg
283 type(ksp_monitor_t) :: ksp_results
284 integer,
optional,
intent(in) :: niter
285 integer :: iter, max_iter
286 real(kind=rp) :: rtr, rnorm
289 if (this%recompute_diag)
then
290 call this%comp_diag(amg, n)
293 if (
present(niter))
then
296 max_iter = this%max_iter
300 associate( w => this%w, r => this%r, d => this%d)
301 do iter = 1, max_iter
304 call amg%matvec(w, x, this%lvl)
312 call add2s2(x, r, this%omega, n)
318 integer,
intent(in) :: lvl
320 character(len=LOG_SIZE) :: log_buf
322 write(log_buf,
'(A8,I2,A28)')
'-- level',lvl,
'-- init smoother: Chebyshev'
323 call neko_log%message(log_buf)
324 write(log_buf,
'(A22,I6)')
'Iterations:',smoo%max_iter
325 call neko_log%message(log_buf)
329 integer,
intent(in) :: lvl
330 real(kind=rp),
intent(in) :: lam
331 character(len=LOG_SIZE) :: log_buf
333 write(log_buf,
'(A12,I2,A29,F12.3)')
'-- AMG level',lvl,
'-- Chebyshev approx. max eig', lam
334 call neko_log%message(log_buf)
Implements the base abstract type for Krylov solvers plus helper types.
type(log_t), public neko_log
Global log stream.
integer, parameter, public log_size
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
integer, parameter, public rp
Global precision used in computations.
Implements smoothers for use with TreeAMG matrix vector product.
subroutine amg_jacobi_init(this, n, lvl, max_iter)
Initialization of Jacobi (this is expensive...)
subroutine amg_cheby_power(this, amg, n)
Power method to approximate largest eigenvalue.
subroutine amg_cheby_init(this, n, lvl, max_iter)
Initialization of chebyshev.
subroutine amg_jacobi_diag(this, amg, n)
SAMPLE MATRIX DIAGONAL VALUES (DO NOT USE, EXPENSIVE)
subroutine amg_jacobi_solve(this, x, f, n, amg, niter)
Jacobi smoother.
subroutine amg_smoo_monitor(lvl, smoo)
subroutine amg_cheby_solve(this, x, f, n, amg, niter)
Chebyshev smoother From Saad's iterative methods textbook.
subroutine amg_cheby_monitor(lvl, lam)
Implements utilities for the TreeAMG hierarchy structure.
Implements the base type for TreeAMG hierarchy structure.
A list of allocatable `bc_t`. Follows the standard interface of lists.
Type for storing initial and final residuals in a Krylov solver.
Type for a TreeAMG hierarchy.
Type for Chebyshev iteration using TreeAMG matvec.
Type for Chebyshev iteration using TreeAMG matvec.