Neko 0.9.99
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tree_amg_smoother.f90
Go to the documentation of this file.
1! Copyright (c) 2024, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use tree_amg
37 use num_types
38 use utils
39 use math
40 use krylov, only : ksp_monitor_t
41 use bc_list, only: bc_list_t
42 use gather_scatter, only : gs_t, gs_op_add
43 use logger, only : neko_log, log_size
44 implicit none
45 private
46
48 type, public :: amg_jacobi_t
49 real(kind=rp), allocatable :: d(:)
50 real(kind=rp), allocatable :: w(:)
51 real(kind=rp), allocatable :: r(:)
52 real(kind=rp) :: omega
53 integer :: lvl
54 integer :: n
55 integer :: max_iter = 10
56 logical :: recompute_diag = .true.
57 contains
58 procedure, pass(this) :: init => amg_jacobi_init
59 procedure, pass(this) :: solve => amg_jacobi_solve
60 procedure, pass(this) :: comp_diag => amg_jacobi_diag
61 end type amg_jacobi_t
62
64 type, public :: amg_cheby_t
65 real(kind=rp), allocatable :: d(:)
66 real(kind=rp), allocatable :: w(:)
67 real(kind=rp), allocatable :: r(:)
68 real(kind=rp) :: tha, dlt
69 integer :: lvl
70 integer :: n
71 integer :: power_its = 250
72 integer :: max_iter = 10
73 logical :: recompute_eigs = .true.
74 contains
75 procedure, pass(this) :: init => amg_cheby_init
76 procedure, pass(this) :: solve => amg_cheby_solve
77 procedure, pass(this) :: comp_eig => amg_cheby_power
78 end type amg_cheby_t
79
80contains
81
86 subroutine amg_cheby_init(this, n, lvl, max_iter)
87 class(amg_cheby_t), intent(inout), target :: this
88 integer, intent(in) :: n
89 integer, intent(in) :: lvl
90 integer, intent(in) :: max_iter
91
92 allocate(this%d(n))
93 allocate(this%w(n))
94 allocate(this%r(n))
95 this%n = n
96 this%lvl = lvl
97 this%max_iter = max_iter
98 this%recompute_eigs = .true.
99
100 call amg_smoo_monitor(lvl,this)
101
102 end subroutine amg_cheby_init
103
104
108 subroutine amg_cheby_power(this, amg, n)
109 class(amg_cheby_t), intent(inout) :: this
110 type(tamg_hierarchy_t), intent(inout) :: amg
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
116 integer :: i
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)
119
120 do i = 1, n
121 !TODO: replace with a better way to initialize power method
122 !call random_number(rn)
123 !d(i) = rn + 10.0_rp
124 d(i) = sin(real(i))
125 end do
126 if (this%lvl .eq. 0) then
127 call gs_h%op(d, n, gs_op_add)!TODO
128 call blst%apply(d, n)
129 end if
130
131 !Power method to get lamba max
132 do i = 1, this%power_its
133 w = 0d0
134 call amg%matvec(w, d, this%lvl)
135
136 if (this%lvl .eq. 0) then
137 !call blst%apply(w, n)
138 wtw = glsc3(w, coef%mult, w, n)
139 else
140 wtw = glsc2(w, w, n)
141 end if
142
143 call cmult2(d, w, 1.0_rp/sqrt(wtw), n)
144
145 if (this%lvl .eq. 0) then
146 !call blst%apply(d, n)
147 end if
148 end do
149
150 w = 0d0
151 call amg%matvec(w, d, this%lvl)
152 if (this%lvl .eq. 0) then
153 !call blst%apply(w, n)
154 end if
155
156 if (this%lvl .eq. 0) then
157 dtw = glsc3(d, coef%mult, w, n)
158 dtd = glsc3(d, coef%mult, d, n)
159 else
160 dtw = glsc2(d, w, n)
161 dtd = glsc2(d, d, n)
162 end if
163 lam = dtw / dtd
164 b = lam * boost
165 a = lam / lam_factor
166 this%tha = (b+a)/2.0_rp
167 this%dlt = (b-a)/2.0_rp
168
169 this%recompute_eigs = .false.
170
171 call amg_cheby_monitor(this%lvl,lam)
172 end associate
173 end subroutine amg_cheby_power
174
181 subroutine amg_cheby_solve(this, x, f, n, amg, niter)
182 class(amg_cheby_t), intent(inout) :: this
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
192
193 if (this%recompute_eigs) then
194 call this%comp_eig(amg, n)
195 end if
196
197 if (present(niter)) then
198 max_iter = niter
199 else
200 max_iter = this%max_iter
201 end if
202
203 associate( w => this%w, r => this%r, d => this%d, blst=>amg%blst)
204 call copy(r, f, n)
205 w = 0d0
206 call amg%matvec(w, x, this%lvl)
207 call sub2(r, w, n)
208
209 thet = this%tha
210 delt = this%dlt
211 s1 = thet / delt
212 rhok = 1.0_rp / s1
213
214 call copy(d, r, n)
215 call cmult(d, 1.0_rp/thet, n)
216
217 do iter = 1, max_iter
218 call add2(x,d,n)
219
220 w = 0d0
221 call amg%matvec(w, d, this%lvl)
222 call sub2(r, w, n)
223
224 rhokp1 = 1.0_rp / (2.0_rp * s1 - rhok)
225 tmp1 = rhokp1 * rhok
226 tmp2 = 2.0_rp * rhokp1 / delt
227 rhok = rhokp1
228
229 call cmult(d, tmp1, n)
230 call add2s2(d, r, tmp2, n)
231 end do
232
233 end associate
234 end subroutine amg_cheby_solve
235
240 subroutine amg_jacobi_init(this, n, lvl, max_iter)
241 class(amg_jacobi_t), intent(inout), target :: this
242 integer, intent(in) :: n
243 integer, intent(in) :: lvl
244 integer, intent(in) :: max_iter
245
246 allocate(this%d(n))
247 allocate(this%w(n))
248 allocate(this%r(n))
249 this%n = n
250 this%lvl = lvl
251 this%max_iter = max_iter
252 this%omega = 0.7_rp
253
254 end subroutine amg_jacobi_init
255
259 subroutine amg_jacobi_diag(this, amg, n)
260 class(amg_jacobi_t), intent(inout) :: this
261 type(tamg_hierarchy_t), intent(inout) :: amg
262 integer, intent(in) :: n
263 real(kind=rp) :: val
264 integer :: i
265 do i = 1, n
266 call tamg_sample_matrix_val(val, amg, this%lvl, i, i)
267 this%d(i) = 1.0_rp / val
268 end do
269 this%recompute_diag = .false.
270 end subroutine amg_jacobi_diag
271
277 subroutine amg_jacobi_solve(this, x, f, n, amg, niter)
278 class(amg_jacobi_t), intent(inout) :: this
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
287 integer :: i
288
289 if (this%recompute_diag) then
290 call this%comp_diag(amg, n)
291 end if
292
293 if (present(niter)) then
294 max_iter = niter
295 else
296 max_iter = this%max_iter
297 end if
298
299 ! x = x + omega * Dinv( f - Ax )
300 associate( w => this%w, r => this%r, d => this%d)
301 do iter = 1, max_iter
302 w = 0.0_rp
304 call amg%matvec(w, x, this%lvl)
306 call copy(r, f, n)
307 call sub2(r, w, n)
308 !print *, iter, "RES", sqrt(glsc2(r,r,n))!>DEBUG
310 call col2(r, d, n)
312 call add2s2(x, r, this%omega, n)
313 end do
314 end associate
315 end subroutine amg_jacobi_solve
316
317 subroutine amg_smoo_monitor(lvl,smoo)
318 integer, intent(in) :: lvl
319 class(amg_cheby_t), intent(in) :: smoo
320 character(len=LOG_SIZE) :: log_buf
321
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)
326 end subroutine amg_smoo_monitor
327
328 subroutine amg_cheby_monitor(lvl,lam)
329 integer, intent(in) :: lvl
330 real(kind=rp), intent(in) :: lam
331 character(len=LOG_SIZE) :: log_buf
332
333 write(log_buf, '(A12,I2,A29,F12.3)') '-- AMG level',lvl,'-- Chebyshev approx. max eig', lam
334 call neko_log%message(log_buf)
335 end subroutine amg_cheby_monitor
336
337end module tree_amg_smoother
double real
Defines a list of bc_t.
Definition bc_list.f90:34
Gather-scatter.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:65
integer, parameter, public log_size
Definition log.f90:42
Definition math.f90:60
subroutine, public cmult2(a, b, c, n)
Multiplication by constant c .
Definition math.f90:700
real(kind=rp) function, public glsc3(a, b, c, n)
Weighted inner product .
Definition math.f90:894
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:875
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
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.
Definition tree_amg.f90:34
Utilities.
Definition utils.f90:35
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:46
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:83
Type for Chebyshev iteration using TreeAMG matvec.
Type for Chebyshev iteration using TreeAMG matvec.