Neko 1.99.3
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, copy
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
63 use logger, only : neko_log, log_size
67 use, intrinsic :: iso_c_binding
68 implicit none
69 private
70
71 type :: tamg_wrk_t
72 integer :: n = -1
73 real(kind=rp), allocatable :: r(:)
74 real(kind=rp), allocatable :: b(:)
75 real(kind=rp), allocatable :: x(:)
76 type(c_ptr) :: r_d = c_null_ptr
77 type(c_ptr) :: b_d = c_null_ptr
78 type(c_ptr) :: x_d = c_null_ptr
79 end type tamg_wrk_t
80
82 type, public :: tamg_solver_t
83 type(tamg_hierarchy_t), allocatable :: amg
84 type(amg_cheby_t), allocatable :: smoo(:)
85 type(tamg_wrk_t), allocatable :: wrk(:)
86 !type(amg_jacobi_t), allocatable :: jsmoo(:)
87 integer :: nlvls
88 integer :: max_iter
89 contains
90 procedure, pass(this) :: init => tamg_mg_init
91 procedure, pass(this) :: solve => tamg_mg_solve
92 procedure, pass(this) :: free => tamg_mg_free
93 procedure, private, pass(this) :: mg_cycle => tamg_mg_cycle
94 procedure, private, pass(this) :: mg_cycle_d => tamg_mg_cycle_d
95 end type tamg_solver_t
96
97contains
98
108 subroutine tamg_mg_init(this, ax, Xh, coef, msh, gs_h, nlvls, blst, &
109 max_iter, cheby_degree)
110 class(tamg_solver_t), intent(inout), target :: this
111 class(ax_t), target, intent(in) :: ax
112 type(space_t), target, intent(in) :: Xh
113 type(coef_t), target, intent(in) :: coef
114 type(mesh_t), target, intent(in) :: msh
115 type(gs_t), target, intent(in) :: gs_h
116 type(bc_list_t), target, intent(in) :: blst
117 integer, intent(in) :: nlvls
118 integer, intent(in) :: max_iter
119 integer, intent(in) :: cheby_degree
120 integer :: lvl, n, mlvl, target_num_aggs
121 integer, allocatable :: agg_nhbr(:,:), nhbr_tmp(:,:)
122 character(len=LOG_SIZE) :: log_buf
123 integer :: glb_min_target_aggs
124 logical :: use_greedy_agg
125
126 call neko_log%section('AMG')
127
128 write(log_buf, '(A28,I2,A8)') 'Creating AMG hierarchy with', &
129 nlvls, 'levels.'
130 call neko_log%message(log_buf)
131
132 allocate( this%amg )
133 call this%amg%init(ax, xh, coef, msh, gs_h, nlvls, blst)
134
135 ! Aggregation
136 use_greedy_agg = .true.
137 ! Create level 1 (neko elements are level 0)
138 call aggregate_finest_level(this%amg, xh%lx, xh%ly, xh%lz, msh%nelv)
139
140 ! Create the remaining levels
141 allocate( agg_nhbr, source = msh%facet_neigh )
142 do mlvl = 2, nlvls-1
143 ! estimate number of aggregates
144 if (use_greedy_agg) then
145 target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 8
146 else
147 target_num_aggs = this%amg%lvl(mlvl-1)%nnodes / 2
148 end if
149
150 glb_min_target_aggs = target_num_aggs
151 call mpi_allreduce(mpi_in_place, glb_min_target_aggs, 1, &
152 mpi_integer, mpi_min, neko_comm)
153 if (glb_min_target_aggs .lt. 4 ) then
154 call neko_warning( &
155 "TAMG: Too many levels. Not enough DOFs for coarsest grid.")
156 this%amg%nlvls = mlvl
157 exit
158 end if
159
160 if (use_greedy_agg) then
161 call print_preagg_info( mlvl, glb_min_target_aggs, 1)
162 call aggregate_greedy( this%amg, mlvl, target_num_aggs, agg_nhbr, nhbr_tmp)
163 else
164 call print_preagg_info( mlvl, glb_min_target_aggs, 2)
165 call aggregate_pairs( this%amg, mlvl, target_num_aggs, agg_nhbr, nhbr_tmp)
166 end if
167
168 agg_nhbr = nhbr_tmp
169 deallocate( nhbr_tmp )
170 end do
171 deallocate( agg_nhbr )
172
173 ! Create the end point
174 call aggregate_end(this%amg, this%amg%nlvls)
175
176 this%max_iter = max_iter
177
178 this%nlvls = this%amg%nlvls
179 if (this%nlvls .gt. this%amg%nlvls) then
180 call neko_error( &
181 "Requested number multigrid levels &
182 & is greater than the initialized AMG levels")
183 end if
184
185 ! Initialize relaxation methods
186 allocate(this%smoo(0:(this%amg%nlvls)))
187 do lvl = 0, this%amg%nlvls-1
188 n = this%amg%lvl(lvl+1)%fine_lvl_dofs
189 call this%smoo(lvl)%init(n, lvl, cheby_degree)
190 end do
191
192 ! Allocate work space on each level
193 allocate(this%wrk(0:(this%amg%nlvls)))
194 do lvl = 0, this%amg%nlvls-1
195 n = this%amg%lvl(lvl+1)%fine_lvl_dofs
196 this%wrk(lvl)%n = n
197 allocate( this%wrk(lvl)%r(n) )
198 allocate( this%wrk(lvl)%b(n) )
199 allocate( this%wrk(lvl)%x(n) )
200 if (neko_bcknd_device .eq. 1) then
201 call device_map( this%wrk(lvl)%r, this%wrk(lvl)%r_d, n)
202 call device_map( this%wrk(lvl)%b, this%wrk(lvl)%b_d, n)
203 call device_map( this%wrk(lvl)%x, this%wrk(lvl)%x_d, n)
204 end if
205 end do
206
207 !allocate(this%jsmoo(0:(this%amg%nlvls)))
208 !do lvl = 0, this%amg%nlvls-1
209 ! n = this%amg%lvl(lvl+1)%fine_lvl_dofs
210 ! call this%jsmoo(lvl)%init(n ,lvl, cheby_degree)
211 !end do
212
213 ! Create index mapping between levels
214 call fill_lvl_map(this%amg)
215
216 call neko_log%end_section()
217
218 end subroutine tamg_mg_init
219
221 subroutine tamg_mg_free(this)
222 class(tamg_solver_t), intent(inout), target :: this
223 integer :: i
224 if (allocated(this%amg)) then
225 call this%amg%free()
226 deallocate(this%amg)
227 end if
228 if (allocated(this%smoo)) then
229 do i = 0, (size(this%smoo)-1)
230 call this%smoo(i)%free()
231 end do
232 deallocate(this%smoo)
233 end if
234 if (allocated(this%wrk)) then
235 do i = 0, (size(this%wrk)-1)
236 if (neko_bcknd_device .eq. 1) then
237 call device_free(this%wrk(i)%r_d)
238 call device_free(this%wrk(i)%b_d)
239 call device_free(this%wrk(i)%x_d)
240 end if
241 if (allocated(this%wrk(i)%r)) deallocate(this%wrk(i)%r)
242 if (allocated(this%wrk(i)%b)) deallocate(this%wrk(i)%b)
243 if (allocated(this%wrk(i)%x)) deallocate(this%wrk(i)%x)
244 end do
245 end if
246 end subroutine tamg_mg_free
247
248
253 subroutine tamg_mg_solve(this, z, r, n)
254 integer, intent(in) :: n
255 class(tamg_solver_t), intent(inout) :: this
256 real(kind=rp), dimension(n), intent(inout) :: z
257 real(kind=rp), dimension(n), intent(inout) :: r
258 type(c_ptr) :: z_d
259 type(c_ptr) :: r_d
260 integer :: iter, max_iter
261 logical :: zero_initial_guess
262
263 max_iter = this%max_iter
264
265 if (neko_bcknd_device .eq. 1) then
266 z_d = device_get_ptr(z)
267 r_d = device_get_ptr(r)
268 ! Zero out the initial guess becuase we do not handle null spaces very well...
269 call device_rzero(this%wrk(0)%x_d, n)
270 call device_copy(this%wrk(0)%b_d, r_d, n)
271 zero_initial_guess = .true.
272 ! Call the amg cycle
273 do iter = 1, max_iter
274 call this%mg_cycle_d(zero_initial_guess)
275 zero_initial_guess = .false.
276 end do
277 call device_copy(z_d, this%wrk(0)%x_d, n)
278 else
279 ! Zero out the initial guess becuase we do not handle null spaces very well...
280 call rzero(this%wrk(0)%x, n)
281 call copy(this%wrk(0)%b, r, n)
282 zero_initial_guess = .true.
283 ! Call the amg cycle
284 do iter = 1, max_iter
285 call this%mg_cycle(zero_initial_guess)
286 zero_initial_guess = .false.
287 end do
288 call copy(z, this%wrk(0)%x, n)
289 end if
290 end subroutine tamg_mg_solve
291
292
296 subroutine tamg_mg_cycle(this, zero_initial_guess)
297 class(tamg_solver_t), intent(inout), target :: this
298 logical, intent(inout) :: zero_initial_guess
299 character(len=2) :: lvl_name
300 integer :: max_lvl, lvl
301
302 max_lvl = this%nlvls-1
303 ! Loop down hierarchy. Fine to coarse
304 do lvl = 0, max_lvl-1
305 write(lvl_name, '(I0)') lvl
306 call profiler_start_region( "AMG_level_" // trim(lvl_name))
307 associate(x => this%wrk(lvl)%x, b => this%wrk(lvl)%b, &
308 r => this%wrk(lvl)%r, n => this%wrk(lvl)%n)
309 !!----------!!
310 !! SMOOTH !!
311 !!----------!!
312 call this%smoo(lvl)%solve(x, b, n, this%amg, &
313 zero_initial_guess)
314 !!----------!!
315 !! Residual !!
316 !!----------!!
317 call calc_resid(r, x, b, this%amg, lvl, n)
318 !!----------!!
319 !! Restrict !!
320 !!----------!!
321 call this%amg%interp_f2c(this%wrk(lvl+1)%b, r, lvl+1)
322
323 call rzero(this%wrk(lvl+1)%x, this%wrk(lvl+1)%n)
324 zero_initial_guess = .true.
325 end associate
326 call profiler_end_region( "AMG_level_" // trim(lvl_name))
327 end do
328 write(lvl_name, '(I0)') max_lvl
329 call profiler_start_region( "AMG_level_" // trim(lvl_name))
330 !!-------------------!!
331 !! Call Coarse solve !!
332 !!-------------------!!
333 call this%smoo(max_lvl)%solve(this%wrk(max_lvl)%x, &
334 this%wrk(max_lvl)%b, this%amg%lvl(max_lvl)%nnodes, this%amg, &
335 zero_initial_guess)
336 call profiler_end_region( "AMG_level_" // trim(lvl_name))
337
338 zero_initial_guess = .false.
339 ! Loop up hierarchy. Coarse to fine
340 do lvl = max_lvl-1, 0, -1
341 write(lvl_name, '(I0)') lvl
342 call profiler_start_region( "AMG_level_" // trim(lvl_name))
343 associate(x => this%wrk(lvl)%x, b => this%wrk(lvl)%b, &
344 r => this%wrk(lvl)%r, n => this%wrk(lvl)%n)
345 !!----------!!
346 !! Project !!
347 !!----------!!
348 call this%amg%interp_c2f(r, this%wrk(lvl+1)%x, lvl+1)
349 !!----------!!
350 !! Correct !!
351 !!----------!!
352 call add2(x, r, n)
353 !!----------!!
354 !! SMOOTH !!
355 !!----------!!
356 call this%smoo(lvl)%solve(x, b, n, this%amg)
357 end associate
358 call profiler_end_region( "AMG_level_" // trim(lvl_name))
359 end do
360 end subroutine tamg_mg_cycle
361
365 subroutine tamg_mg_cycle_d(this, zero_initial_guess)
366 class(tamg_solver_t), intent(inout), target :: this
367 logical, intent(inout) :: zero_initial_guess
368 character(len=2) :: lvl_name
369 integer :: max_lvl, lvl
370
371 max_lvl = this%nlvls-1
372 ! Loop down hierarchy. Fine to coarse
373 do lvl = 0, max_lvl-1
374 write(lvl_name, '(I0)') lvl
375 call profiler_start_region( "AMG_level_" // trim(lvl_name))
376 associate(x => this%wrk(lvl)%x, x_d => this%wrk(lvl)%x_d, &
377 b => this%wrk(lvl)%b, b_d => this%wrk(lvl)%b_d, &
378 r => this%wrk(lvl)%r, r_d => this%wrk(lvl)%r_d, &
379 n => this%wrk(lvl)%n)
380 !!----------!!
381 !! SMOOTH !!
382 !!----------!!
383 call this%smoo(lvl)%device_solve(x, b, x_d, b_d, n, this%amg, &
384 zero_initial_guess)
385 !!----------!!
386 !! Residual !!
387 !!----------!!
388 call this%amg%device_matvec(r, x, r_d, x_d, lvl)
389 call device_sub3(r_d, b_d, r_d, n)
390 !!----------!!
391 !! Restrict !!
392 !!----------!!
393 call this%amg%interp_f2c_d(this%wrk(lvl+1)%b_d, r_d, lvl+1)
394
395 call device_rzero(this%wrk(lvl+1)%x_d, this%wrk(lvl+1)%n)
396 zero_initial_guess = .true.
397 end associate
398 call profiler_end_region( "AMG_level_" // trim(lvl_name))
399 end do
400 write(lvl_name, '(I0)') max_lvl
401 call profiler_start_region( "AMG_level_" // trim(lvl_name))
402 !!-------------------!!
403 !! Call Coarse solve !!
404 !!-------------------!!
405 call this%smoo(max_lvl)%device_solve( &
406 this%wrk(max_lvl)%x, this%wrk(max_lvl)%b, &
407 this%wrk(max_lvl)%x_d, this%wrk(max_lvl)%b_d, &
408 this%amg%lvl(max_lvl)%nnodes, this%amg, &
409 zero_initial_guess)
410 call profiler_end_region( "AMG_level_" // trim(lvl_name))
411
412 zero_initial_guess = .false.
413 ! Loop up hierarchy. Coarse to fine
414 do lvl = max_lvl-1, 0, -1
415 write(lvl_name, '(I0)') lvl
416 call profiler_start_region( "AMG_level_" // trim(lvl_name))
417 associate(x => this%wrk(lvl)%x, x_d => this%wrk(lvl)%x_d, &
418 b => this%wrk(lvl)%b, b_d => this%wrk(lvl)%b_d, &
419 r => this%wrk(lvl)%r, r_d => this%wrk(lvl)%r_d, &
420 n => this%wrk(lvl)%n)
421 !!----------!!
422 !! Project !!
423 !!----------!!
424 call this%amg%interp_c2f_d(r_d, this%wrk(lvl+1)%x_d, lvl+1, r)
425 !!----------!!
426 !! Correct !!
427 !!----------!!
428 call device_add2(x_d, r_d, n)
429 !!----------!!
430 !! SMOOTH !!
431 !!----------!!
432 call this%smoo(lvl)%device_solve(x, b, x_d, b_d, n, this%amg)
433 end associate
434 call profiler_end_region( "AMG_level_" // trim(lvl_name))
435 end do
436 end subroutine tamg_mg_cycle_d
437
438
446 subroutine calc_resid(r, x, b, amg, lvl, n)
447 integer, intent(in) :: n
448 real(kind=rp), intent(inout) :: r(n)
449 real(kind=rp), intent(inout) :: x(n)
450 real(kind=rp), intent(inout) :: b(n)
451 type(tamg_hierarchy_t), intent(inout) :: amg
452 integer, intent(in) :: lvl
453 integer :: i
454 call amg%matvec(r, x, lvl)
455 call sub3(r, b, r, n)
456 end subroutine calc_resid
457
458
459 subroutine print_preagg_info(lvl, nagg, agg_type)
460 integer, intent(in) :: lvl, nagg, agg_type
461 character(len=LOG_SIZE) :: log_buf
462 !TODO: calculate min and max agg size
463 if (agg_type .eq. 1) then
464 write(log_buf, '(A8,I2,A31)') '-- level', lvl, &
465 '-- Calling Greedy Aggregation'
466 else if (agg_type .eq. 2) then
467 write(log_buf, '(A8,I2,A33)') '-- level', lvl, &
468 '-- Calling Pairwise Aggregation'
469 else
470 write(log_buf, '(A8,I2,A31)') '-- level', lvl, &
471 '-- UNKNOWN Aggregation'
472 end if
473 call neko_log%message(log_buf)
474 write(log_buf, '(A33,I6)') 'Target Aggregates:', nagg
475 call neko_log%message(log_buf)
476 end subroutine print_preagg_info
477
478 subroutine print_resid_info(r, x, b, r_d, x_d, b_d, amg, lvl, n)
479 integer, intent(in) :: lvl, n
480 real(kind=rp), intent(inout) :: r(n)
481 real(kind=rp), intent(inout) :: x(n)
482 real(kind=rp), intent(inout) :: b(n)
483 type(c_ptr) :: r_d
484 type(c_ptr) :: x_d
485 type(c_ptr) :: b_d
486 type(tamg_hierarchy_t), intent(inout) :: amg
487 real(kind=rp) :: val
488 character(len=LOG_SIZE) :: log_buf
489
490 call amg%device_matvec(r, x, r_d, x_d, lvl)
491 call device_sub3(r_d, b_d, r_d, n)
492 val = device_glsc2(r_d, r_d, n)
493
494 write(log_buf, '(A33,I6,F12.6)') 'tAMG resid:', lvl, val
495 call neko_log%message(log_buf)
496 end subroutine print_resid_info
497
500 subroutine fill_lvl_map(amg)
501 type(tamg_hierarchy_t), intent(inout) :: amg
502 integer :: i, j, k, l, nid, n
503 do j = 1, amg%lvl(1)%nnodes
504 do k = 1, amg%lvl(1)%nodes(j)%ndofs
505 nid = amg%lvl(1)%nodes(j)%dofs(k)
506 amg%lvl(1)%map_finest2lvl(nid) = amg%lvl(1)%nodes(j)%gid
507 end do
508 end do
509 n = size(amg%lvl(1)%map_finest2lvl)
510 do l = 2, amg%nlvls
511 do i = 1, n
512 nid = amg%lvl(l-1)%map_finest2lvl(i)
513 do j = 1, amg%lvl(l)%nnodes
514 do k = 1, amg%lvl(l)%nodes(j)%ndofs
515 if (nid .eq. amg%lvl(l)%nodes(j)%dofs(k)) then
516 amg%lvl(l)%map_finest2lvl(i) = amg%lvl(l)%nodes(j)%gid
517 end if
518 end do
519 end do
520 end do
521 end do
522 if (neko_bcknd_device .eq. 1) then
523 do l = 1, amg%nlvls
524 amg%lvl(l)%map_finest2lvl(0) = n
525 call device_memcpy( amg%lvl(l)%map_finest2lvl, &
526 amg%lvl(l)%map_finest2lvl_d, n, &
527 host_to_device, .true.)
528 call device_memcpy( amg%lvl(l)%map_f2c, &
529 amg%lvl(l)%map_f2c_d, amg%lvl(l)%fine_lvl_dofs+1, &
530 host_to_device, .true.)
531 end do
532 end if
533 end subroutine fill_lvl_map
534end 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:107
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_copy(a_d, b_d, n, strm)
Copy a 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:225
Gather-scatter.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
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 copy(a, b, n)
Copy a vector .
Definition math.f90:249
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
Profiling interface.
Definition profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:107
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:
subroutine calc_resid(r, x, b, amg, lvl, n)
Wrapper function to calculate residyal.
subroutine tamg_mg_cycle_d(this, zero_initial_guess)
multigrid cycle for the TreeAMG solver object on device
subroutine tamg_mg_free(this)
free tree amg solver object
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 tamg_mg_cycle(this, zero_initial_guess)
multigrid cycle 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:255
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
Definition tree_amg.f90:186
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:61
The function space for the SEM solution fields.
Definition space.f90:63
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:87
Type for the TreeAMG solver.
Type for Chebyshev iteration using TreeAMG matvec.