Neko 1.99.1
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tree_amg_multigrid.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!
45 use num_types, only: rp
46 use utils, only : neko_error, neko_warning
47 use math, only : add2, rzero, glsc2, sub3, col2
50 use comm
51 use mpi_f08, only: mpi_allreduce, mpi_min, mpi_in_place, mpi_integer
52 use coefs, only : coef_t
53 use mesh, only : mesh_t
54 use space, only : space_t
55 use ax_product, only: ax_t
56 use bc_list, only : bc_list_t
57 use gather_scatter, only : gs_t, gs_op_add
62 use logger, only : neko_log, log_size
66 use, intrinsic :: iso_c_binding
67 implicit none
68 private
69
70 type :: tamg_wrk_t
71 real(kind=rp), allocatable :: r(:)
72 real(kind=rp), allocatable :: rc(:)
73 real(kind=rp), allocatable :: tmp(:)
74 type(c_ptr) :: r_d = c_null_ptr
75 type(c_ptr) :: rc_d = c_null_ptr
76 type(c_ptr) :: tmp_d = c_null_ptr
77 end type tamg_wrk_t
78
80 type, public :: tamg_solver_t
81 type(tamg_hierarchy_t), allocatable :: amg
82 type(amg_cheby_t), allocatable :: smoo(:)
83 type(tamg_wrk_t), allocatable :: wrk(:)
84 !type(amg_jacobi_t), allocatable :: jsmoo(:)
85 integer :: nlvls
86 integer :: max_iter
87 contains
88 procedure, pass(this) :: init => tamg_mg_init
89 procedure, pass(this) :: solve => tamg_mg_solve
90 end type tamg_solver_t
91
92contains
93
103 subroutine tamg_mg_init(this, ax, Xh, coef, msh, gs_h, nlvls, blst, &
104 max_iter, cheby_degree)
105 class(tamg_solver_t), intent(inout), target :: this
106 class(ax_t), target, intent(in) :: ax
107 type(space_t), target, intent(in) :: Xh
108 type(coef_t), target, intent(in) :: coef
109 type(mesh_t), target, intent(in) :: msh
110 type(gs_t), target, intent(in) :: gs_h
111 type(bc_list_t), target, intent(in) :: blst
112 integer, intent(in) :: nlvls
113 integer, intent(in) :: max_iter
114 integer, intent(in) :: cheby_degree
115 integer :: lvl, n, mlvl, target_num_aggs
116 integer, allocatable :: agg_nhbr(:,:), nhbr_tmp(:,:)
117 character(len=LOG_SIZE) :: log_buf
118 integer :: glb_min_target_aggs
119 logical :: use_greedy_agg
120
121 call neko_log%section('AMG')
122
123 write(log_buf, '(A28,I2,A8)') 'Creating AMG hierarchy with', &
124 nlvls, 'levels.'
125 call neko_log%message(log_buf)
126
127 allocate( this%amg )
128 call this%amg%init(ax, xh, coef, msh, gs_h, nlvls, blst)
129
130 ! Aggregation
131 use_greedy_agg = .true.
132 ! Create level 1 (neko elements are level 0)
133 call aggregate_finest_level(this%amg, xh%lx, xh%ly, xh%lz, msh%nelv)
134
135 ! Create the remaining levels
136 allocate( agg_nhbr, source = msh%facet_neigh )
137 do mlvl = 2, nlvls-1
138 ! estimate number of aggregates
139 if (use_greedy_agg) then
140 target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 8
141 else
142 target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 2
143 end if
144
145 glb_min_target_aggs = target_num_aggs
146 call mpi_allreduce(mpi_in_place, glb_min_target_aggs, 1, &
147 mpi_integer, mpi_min, neko_comm)
148 if (glb_min_target_aggs .lt. 4 ) then
149 call neko_warning( &
150 "TAMG: Too many levels. Not enough DOFs for coarsest grid.")
151 this%amg%nlvls = mlvl
152 exit
153 end if
154
155 if (use_greedy_agg) then
156 call print_preagg_info( mlvl, glb_min_target_aggs, 1)
157 call aggregate_greedy( this%amg, mlvl, target_num_aggs, agg_nhbr, nhbr_tmp)
158 else
159 call print_preagg_info( mlvl, glb_min_target_aggs, 2)
160 call aggregate_pairs( this%amg, mlvl, target_num_aggs, agg_nhbr, nhbr_tmp)
161 end if
162
163 agg_nhbr = nhbr_tmp
164 deallocate( nhbr_tmp )
165 end do
166 deallocate( agg_nhbr )
167
168 ! Create the end point
169 call aggregate_end(this%amg, this%amg%nlvls)
170
171 this%max_iter = max_iter
172
173 this%nlvls = this%amg%nlvls
174 if (this%nlvls .gt. this%amg%nlvls) then
175 call neko_error( &
176 "Requested number multigrid levels &
177 & is greater than the initialized AMG levels")
178 end if
179
180 ! Initialize relaxation methods
181 allocate(this%smoo(0:(this%amg%nlvls)))
182 do lvl = 0, this%amg%nlvls-1
183 n = this%amg%lvl(lvl+1)%fine_lvl_dofs
184 call this%smoo(lvl)%init(n, lvl, cheby_degree)
185 end do
186
187 ! Allocate work space on each level
188 allocate(this%wrk(this%amg%nlvls))
189 do lvl = 1, this%amg%nlvls
190 n = this%amg%lvl(lvl)%fine_lvl_dofs
191 allocate( this%wrk(lvl)%r(n) )
192 allocate( this%wrk(lvl)%rc(n) )
193 allocate( this%wrk(lvl)%tmp(n) )
194 if (neko_bcknd_device .eq. 1) then
195 call device_map( this%wrk(lvl)%r, this%wrk(lvl)%r_d, n)
196 call device_map( this%wrk(lvl)%rc, this%wrk(lvl)%rc_d, n)
197 call device_map( this%wrk(lvl)%tmp, this%wrk(lvl)%tmp_d, n)
198 end if
199 end do
200
201 !allocate(this%jsmoo(0:(this%amg%nlvls)))
202 !do lvl = 0, this%amg%nlvls-1
203 ! n = this%amg%lvl(lvl+1)%fine_lvl_dofs
204 ! call this%jsmoo(lvl)%init(n ,lvl, cheby_degree)
205 !end do
206
207 ! Create index mapping between levels
208 call fill_lvl_map(this%amg)
209
210 call neko_log%end_section()
211
212 end subroutine tamg_mg_init
213
214
219 subroutine tamg_mg_solve(this, z, r, n)
220 integer, intent(in) :: n
221 class(tamg_solver_t), intent(inout) :: this
222 real(kind=rp), dimension(n), intent(inout) :: z
223 real(kind=rp), dimension(n), intent(inout) :: r
224 type(c_ptr) :: z_d
225 type(c_ptr) :: r_d
226 integer :: iter, max_iter
227 logical :: zero_initial_guess
228
229 max_iter = this%max_iter
230
231 if (neko_bcknd_device .eq. 1) then
232 z_d = device_get_ptr(z)
233 r_d = device_get_ptr(r)
234 ! Zero out the initial guess becuase we do not handle null spaces very well...
235 call device_rzero(z_d, n)
236 zero_initial_guess = .true.
237 ! Call the amg cycle
238 do iter = 1, max_iter
239 call tamg_mg_cycle_d(z, r, z_d, r_d, n, 0, this%amg, this, &
240 zero_initial_guess)
241 zero_initial_guess = .false.
242 end do
243 else
244 ! Zero out the initial guess becuase we do not handle null spaces very well...
245 call rzero(z, n)
246 zero_initial_guess = .true.
247 ! Call the amg cycle
248 do iter = 1, max_iter
249 call tamg_mg_cycle(z, r, n, 0, this%amg, this, &
250 zero_initial_guess)
251 zero_initial_guess = .false.
252 end do
253 end if
254 end subroutine tamg_mg_solve
255
256
264 recursive subroutine tamg_mg_cycle(x, b, n, lvl, amg, mgstuff, &
265 zero_initial_guess)
266 integer, intent(in) :: n
267 real(kind=rp), intent(inout) :: x(n)
268 real(kind=rp), intent(inout) :: b(n)
269 type(tamg_hierarchy_t), intent(inout) :: amg
270 type(tamg_solver_t), intent(inout) :: mgstuff
271 logical, intent(in) :: zero_initial_guess
272 integer, intent(in) :: lvl
273 real(kind=rp) :: r(n)
274 real(kind=rp) :: rc(n)
275 real(kind=rp) :: tmp(n)
276 integer :: iter, num_iter
277 integer :: max_lvl
278 integer :: i, cyt
279 max_lvl = mgstuff%nlvls-1
280 !!----------!!
281 !! SMOOTH !!
282 !!----------!!
283 call mgstuff%smoo(lvl)%solve(x, b, n, amg, &
284 zero_initial_guess)
285 if (lvl .eq. max_lvl) then
286 return
287 end if
288 !!----------!!
289 !! Residual !!
290 !!----------!!
291 call calc_resid(r, x, b, amg, lvl, n)
292 !!----------!!
293 !! Restrict !!
294 !!----------!!
295 call amg%interp_f2c(rc, r, lvl+1)
296 !!-------------------!!
297 !! Call Coarse solve !!
298 !!-------------------!!
299 call rzero(tmp, n)
300 call tamg_mg_cycle(tmp, rc, amg%lvl(lvl+1)%nnodes, lvl+1, amg, mgstuff, &
301 .true.)
302 !!----------!!
303 !! Project !!
304 !!----------!!
305 call amg%interp_c2f(r, tmp, lvl+1)
306 !!----------!!
307 !! Correct !!
308 !!----------!!
309 call add2(x, r, n)
310 !!----------!!
311 !! SMOOTH !!
312 !!----------!!
313 call mgstuff%smoo(lvl)%solve(x,b, n, amg)
314 end subroutine tamg_mg_cycle
315
323 recursive subroutine tamg_mg_cycle_d(x, b, x_d, b_d, n, lvl, amg, mgstuff, &
324 zero_initial_guess)
325 integer, intent(in) :: n
326 real(kind=rp), intent(inout) :: x(n)
327 real(kind=rp), intent(inout) :: b(n)
328 type(c_ptr) :: x_d
329 type(c_ptr) :: b_d
330 type(tamg_hierarchy_t), intent(inout) :: amg
331 type(tamg_solver_t), intent(inout) :: mgstuff
332 logical, intent(in) :: zero_initial_guess
333 integer, intent(in) :: lvl
334 integer :: iter, num_iter
335 integer :: max_lvl
336 integer :: i, cyt
337 max_lvl = mgstuff%nlvls-1
338 !!----------!!
339 !! SMOOTH !!
340 !!----------!!
341 call mgstuff%smoo(lvl)%device_solve(x, b, x_d, b_d, n, amg, &
342 zero_initial_guess)
343 if (lvl .eq. max_lvl) then
344 return
345 end if
346 associate( r => mgstuff%wrk(lvl+1)%r, r_d => mgstuff%wrk(lvl+1)%r_d, &
347 rc => mgstuff%wrk(lvl+1)%rc, rc_d => mgstuff%wrk(lvl+1)%rc_d, &
348 tmp => mgstuff%wrk(lvl+1)%tmp, tmp_d => mgstuff%wrk(lvl+1)%tmp_d )
349 !!----------!!
350 !! Residual !!
351 !!----------!!
352 call amg%device_matvec(r, x, r_d, x_d, lvl)
353 call device_sub3(r_d, b_d, r_d, n)
354 !!----------!!
355 !! Restrict !!
356 !!----------!!
357 call amg%interp_f2c_d(rc_d, r_d, lvl+1)
358 !!-------------------!!
359 !! Call Coarse solve !!
360 !!-------------------!!
361 call device_rzero(tmp_d, n)
362 call tamg_mg_cycle_d(tmp, rc, tmp_d, rc_d, &
363 amg%lvl(lvl+1)%nnodes, lvl+1, amg, mgstuff, .true.)
364 !!----------!!
365 !! Project !!
366 !!----------!!
367 call amg%interp_c2f_d(r_d, tmp_d, lvl+1, r)
368 !!----------!!
369 !! Correct !!
370 !!----------!!
371 call device_add2(x_d, r_d, n)
372 !!----------!!
373 !! SMOOTH !!
374 !!----------!!
375 call mgstuff%smoo(lvl)%device_solve(x, b, x_d, b_d, n, amg)
376 end associate
377 end subroutine tamg_mg_cycle_d
378
379
387 subroutine calc_resid(r, x, b, amg, lvl, n)
388 integer, intent(in) :: n
389 real(kind=rp), intent(inout) :: r(n)
390 real(kind=rp), intent(inout) :: x(n)
391 real(kind=rp), intent(inout) :: b(n)
392 type(tamg_hierarchy_t), intent(inout) :: amg
393 integer, intent(in) :: lvl
394 integer :: i
395 call amg%matvec(r, x, lvl)
396 call sub3(r, b, r, n)
397 end subroutine calc_resid
398
399
400 subroutine print_preagg_info(lvl, nagg, agg_type)
401 integer, intent(in) :: lvl, nagg, agg_type
402 character(len=LOG_SIZE) :: log_buf
403 !TODO: calculate min and max agg size
404 if (agg_type .eq. 1) then
405 write(log_buf, '(A8,I2,A31)') '-- level', lvl, &
406 '-- Calling Greedy Aggregation'
407 else if (agg_type .eq. 2) then
408 write(log_buf, '(A8,I2,A33)') '-- level', lvl, &
409 '-- Calling Pairwise Aggregation'
410 else
411 write(log_buf, '(A8,I2,A31)') '-- level', lvl, &
412 '-- UNKNOWN Aggregation'
413 end if
414 call neko_log%message(log_buf)
415 write(log_buf, '(A33,I6)') 'Target Aggregates:', nagg
416 call neko_log%message(log_buf)
417 end subroutine print_preagg_info
418
419 subroutine print_resid_info(r, x, b, r_d, x_d, b_d, amg, lvl, n)
420 integer, intent(in) :: lvl, n
421 real(kind=rp), intent(inout) :: r(n)
422 real(kind=rp), intent(inout) :: x(n)
423 real(kind=rp), intent(inout) :: b(n)
424 type(c_ptr) :: r_d
425 type(c_ptr) :: x_d
426 type(c_ptr) :: b_d
427 type(tamg_hierarchy_t), intent(inout) :: amg
428 real(kind=rp) :: val
429 character(len=LOG_SIZE) :: log_buf
430
431 call amg%device_matvec(r, x, r_d, x_d, lvl)
432 call device_sub3(r_d, b_d, r_d, n)
433 val = device_glsc2(r_d, r_d, n)
434
435 write(log_buf, '(A33,I6,F12.6)') 'tAMG resid:', lvl, val
436 call neko_log%message(log_buf)
437 end subroutine print_resid_info
438
441 subroutine fill_lvl_map(amg)
442 type(tamg_hierarchy_t), intent(inout) :: amg
443 integer :: i, j, k, l, nid, n
444 do j = 1, amg%lvl(1)%nnodes
445 do k = 1, amg%lvl(1)%nodes(j)%ndofs
446 nid = amg%lvl(1)%nodes(j)%dofs(k)
447 amg%lvl(1)%map_finest2lvl(nid) = amg%lvl(1)%nodes(j)%gid
448 end do
449 end do
450 n = size(amg%lvl(1)%map_finest2lvl)
451 do l = 2, amg%nlvls
452 do i = 1, n
453 nid = amg%lvl(l-1)%map_finest2lvl(i)
454 do j = 1, amg%lvl(l)%nnodes
455 do k = 1, amg%lvl(l)%nodes(j)%ndofs
456 if (nid .eq. amg%lvl(l)%nodes(j)%dofs(k)) then
457 amg%lvl(l)%map_finest2lvl(i) = amg%lvl(l)%nodes(j)%gid
458 end if
459 end do
460 end do
461 end do
462 end do
463 if (neko_bcknd_device .eq. 1) then
464 do l = 1, amg%nlvls
465 amg%lvl(l)%map_finest2lvl(0) = n
466 call device_memcpy( amg%lvl(l)%map_finest2lvl, &
467 amg%lvl(l)%map_finest2lvl_d, n, &
468 host_to_device, .true.)
469 call device_memcpy( amg%lvl(l)%map_f2c, &
470 amg%lvl(l)%map_f2c_d, amg%lvl(l)%fine_lvl_dofs+1, &
471 host_to_device, .true.)
472 end do
473 end if
474 end subroutine fill_lvl_map
475end module tree_amg_multigrid
__device__ T solve(const T u, const T y, const T guess, const T nu, const T kappa, const T B)
Return the device pointer for an associated Fortran array.
Definition device.F90:101
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
Copy data between host and device (or device and device)
Definition device.F90:71
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
subroutine, public device_sub3(a_d, b_d, c_d, n, strm)
Vector subtraction .
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition device.F90:219
Gather-scatter.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:70
integer, parameter, public log_size
Definition log.f90:46
Definition math.f90:60
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public sub3(a, b, c, n)
Vector subtraction .
Definition math.f90:781
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:205
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Defines a function space.
Definition space.f90:34
Implements an aggregation for TreeAMG hierarchy structure.
subroutine aggregate_end(tamg, lvl_id)
Aggregate all dofs to a single point to form a tree-like structure.
subroutine aggregate_pairs(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
Aggregates pairs of dofs based on adjacent dofs.
subroutine aggregate_greedy(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
Aggregates dofs based on adjacent dofs.
subroutine aggregate_finest_level(tamg, lx, ly, lz, ne)
Aggregaiton on finest level Aggregates all dofs in an element into a single aggregate.
Implements multigrid using the TreeAMG hierarchy structure. USE:
recursive subroutine tamg_mg_cycle(x, b, n, lvl, amg, mgstuff, zero_initial_guess)
Recrsive multigrid cycle for the TreeAMG solver object.
subroutine calc_resid(r, x, b, amg, lvl, n)
Wrapper function to calculate residyal.
recursive subroutine tamg_mg_cycle_d(x, b, x_d, b_d, n, lvl, amg, mgstuff, zero_initial_guess)
Recrsive multigrid cycle for the TreeAMG solver object on device.
subroutine fill_lvl_map(amg)
Create index mapping between levels and directly to finest level.
subroutine print_resid_info(r, x, b, r_d, x_d, b_d, amg, lvl, n)
subroutine tamg_mg_solve(this, z, r, n)
Solver function for the TreeAMG solver object.
subroutine print_preagg_info(lvl, nagg, agg_type)
subroutine tamg_mg_init(this, ax, xh, coef, msh, gs_h, nlvls, blst, max_iter, cheby_degree)
Initialization of the TreeAMG multigrid solver.
Implements smoothers for use with TreeAMG matrix vector product.
Implements the base type for TreeAMG hierarchy structure.
Definition tree_amg.f90:34
subroutine, public tamg_node_init(node, gid, ndofs)
Initialization of a TreeAMG tree node.
Definition tree_amg.f90:203
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
Definition tree_amg.f90:171
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:55
The function space for the SEM solution fields.
Definition space.f90:63
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:94
Type for the TreeAMG solver.
Type for Chebyshev iteration using TreeAMG matvec.