Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
tree_amg.f90
Go to the documentation of this file.
1! Copyright (c) 2024-2026, 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 num_types, only : rp
36 use utils, only : neko_error
37 use math, only : rzero, col2
40 use coefs, only : coef_t
41 use mesh, only : mesh_t
42 use space, only : space_t
43 use ax_product, only: ax_t
44 use bc_list, only: bc_list_t
45 use gather_scatter, only : gs_t, gs_op_add
46 use device, only: device_map, device_unmap, &
49 use, intrinsic :: iso_c_binding
50 implicit none
51 private
52
54 type, private :: tamg_node_t
55 logical :: isleaf = .true.
56 integer :: gid = -1
57 integer :: lvl = -1
58 integer :: ndofs = 0
59 integer, allocatable :: dofs(:)
60 real(kind=rp) :: xyz(3)
61 real(kind=rp), allocatable :: interp_r(:)
62 real(kind=rp), allocatable :: interp_p(:)
63 contains
64 procedure, pass(this) :: free => node_free
65 end type tamg_node_t
66
68 type, private :: tamg_lvl_t
69 integer :: lvl = -1
70 integer :: nnodes = 0
71 type(tamg_node_t), allocatable :: nodes(:)
72 integer :: fine_lvl_dofs = 0
73 real(kind=rp), allocatable :: wrk_in(:)
74 type(c_ptr) :: wrk_in_d = c_null_ptr
75 real(kind=rp), allocatable :: wrk_out(:)
76 type(c_ptr) :: wrk_out_d = c_null_ptr
77 integer, allocatable :: map_finest2lvl(:)
78 type(c_ptr) :: map_finest2lvl_d = c_null_ptr
79 !--!
80 integer, allocatable :: map_f2c(:)
81 type(c_ptr) :: map_f2c_d = c_null_ptr
82 contains
83 procedure, pass(this) :: free => lvl_free
84 end type tamg_lvl_t
85
87 type, public :: tamg_hierarchy_t
89 integer :: nlvls
91 type(tamg_lvl_t), allocatable :: lvl(:)
92
94 class(ax_t), pointer :: ax
95 type(mesh_t), pointer :: msh
96 type(space_t), pointer :: xh
97 type(coef_t), pointer :: coef
98 type(gs_t), pointer :: gs_h
99 type(bc_list_t), pointer :: blst
100
101 contains
102 procedure, pass(this) :: init => tamg_init
103 procedure, pass(this) :: free => tamg_free
104 procedure, pass(this) :: matvec => tamg_matvec
105 procedure, pass(this) :: matvec_impl => tamg_matvec_impl
106 procedure, pass(this) :: interp_f2c => tamg_restriction_operator
107 procedure, pass(this) :: interp_c2f => tamg_prolongation_operator
108 procedure, pass(this) :: interp_f2c_d => tamg_device_restriction_operator
109 procedure, pass(this) :: interp_c2f_d => tamg_device_prolongation_operator
110 procedure, pass(this) :: device_matvec => tamg_device_matvec_flat_impl
111 end type tamg_hierarchy_t
112
114
115contains
116
125 subroutine tamg_init(this, ax, Xh, coef, msh, gs_h, nlvls, blst)
126 class(tamg_hierarchy_t), target, intent(inout) :: this
127 class(ax_t), target, intent(in) :: ax
128 type(space_t),target, intent(in) :: Xh
129 type(coef_t), target, intent(in) :: coef
130 type(mesh_t), target, intent(in) :: msh
131 type(gs_t), target, intent(in) :: gs_h
132 integer, intent(in) :: nlvls
133 type(bc_list_t), target, intent(in) :: blst
134 integer :: i, n
135
136 this%ax => ax
137 this%msh => msh
138 this%Xh => xh
139 this%coef => coef
140 this%gs_h => gs_h
141 this%blst => blst
142
143 if (nlvls .lt. 2) then
144 call neko_error("Need to request at least two multigrid levels.")
145 end if
146
147 this%nlvls = nlvls
148 allocate( this%lvl(this%nlvls) )
149
150 do i = 1, nlvls
151 allocate( this%lvl(i)%map_finest2lvl( 0:coef%dof%size() ))
152 if (neko_bcknd_device .eq. 1) then
153 call device_map(this%lvl(i)%map_finest2lvl, this%lvl(i)%map_finest2lvl_d, coef%dof%size()+1)
154 end if
155 end do
156
157 end subroutine tamg_init
158
160 subroutine tamg_free(this)
161 class(tamg_hierarchy_t), intent(inout) :: this
162 integer :: i
163 if (allocated(this%lvl)) then
164 ! using size() instead of this%nlvls since
165 ! this%nlvls may be less than the allocated number of levels due to early
166 ! termination of aggregation
167 do i = 1, size(this%lvl)
168 call this%lvl(i)%free()
169 end do
170 deallocate(this%lvl)
171 end if
172 nullify(this%ax)
173 nullify(this%msh)
174 nullify(this%Xh)
175 nullify(this%coef)
176 nullify(this%gs_h)
177 nullify(this%blst)
178 end subroutine tamg_free
179
185 subroutine tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
186 type(tamg_lvl_t), intent(inout) :: tamg_lvl
187 integer, intent(in) :: lvl
188 integer, intent(in) :: nnodes
189 integer, intent(in) :: ndofs
190
191 tamg_lvl%lvl = lvl
192 tamg_lvl%nnodes = nnodes
193 allocate( tamg_lvl%nodes(tamg_lvl%nnodes) )
194 allocate( tamg_lvl%map_f2c(0:ndofs) )
195 if (neko_bcknd_device .eq. 1) then
196 call device_map(tamg_lvl%map_f2c, tamg_lvl%map_f2c_d, ndofs+1)
197 end if
198
199 tamg_lvl%fine_lvl_dofs = ndofs
200 allocate( tamg_lvl%wrk_in( ndofs ) )
201 allocate( tamg_lvl%wrk_out( ndofs ) )
202 if (neko_bcknd_device .eq. 1) then
203 call device_map(tamg_lvl%wrk_in, tamg_lvl%wrk_in_d, ndofs)
204 call device_cfill(tamg_lvl%wrk_in_d, 0.0_rp, ndofs)
205 call device_map(tamg_lvl%wrk_out, tamg_lvl%wrk_out_d, ndofs)
206 call device_cfill(tamg_lvl%wrk_out_d, 0.0_rp, ndofs)
207 end if
208 end subroutine tamg_lvl_init
209
211 subroutine lvl_free(this)
212 class(tamg_lvl_t), intent(inout) :: this
213 integer :: i
214
215 if (allocated(this%nodes)) then
216 do i = 1, this%nnodes
217 call this%nodes(i)%free()
218 end do
219 deallocate(this%nodes)
220 end if
221 if (allocated(this%wrk_in)) then
222 if (neko_bcknd_device .eq. 1 .and. c_associated(this%wrk_in_d)) then
223 call device_unmap(this%wrk_in, this%wrk_in_d)
224 end if
225 deallocate(this%wrk_in)
226 end if
227 if (allocated(this%wrk_out)) then
228 if (neko_bcknd_device .eq. 1 .and. c_associated(this%wrk_out_d)) then
229 call device_unmap(this%wrk_out, this%wrk_out_d)
230 end if
231 deallocate(this%wrk_out)
232 end if
233 if (allocated(this%map_f2c)) then
234 if (neko_bcknd_device .eq. 1 .and. c_associated(this%map_f2c_d)) then
235 call device_unmap(this%map_f2c, this%map_f2c_d)
236 end if
237 deallocate(this%map_f2c)
238 end if
239 if (allocated(this%map_finest2lvl)) then
240 if (neko_bcknd_device .eq. 1 .and. c_associated(this%map_finest2lvl_d)) then
241 call device_unmap(this%map_finest2lvl, this%map_finest2lvl_d)
242 end if
243 deallocate(this%map_finest2lvl)
244 end if
245 this%nnodes = 0
246 this%lvl = -1
247 this%fine_lvl_dofs = 0
248 end subroutine lvl_free
249
254 subroutine tamg_node_init(node, gid, ndofs)
255 type(tamg_node_t), intent(inout) :: node
256 integer, intent(in) :: gid
257 integer, intent(in) :: ndofs
258
259 node%gid = gid
260 node%ndofs = ndofs
261 allocate( node%dofs( node%ndofs) )
262 node%dofs = -1
263 allocate( node%interp_r( node%ndofs) )
264 allocate( node%interp_p( node%ndofs) )
265 node%interp_r = 1.0_rp
266 node%interp_p = 1.0_rp
267 end subroutine tamg_node_init
268
270 subroutine node_free(this)
271 class(tamg_node_t), intent(inout) :: this
272
273 if (allocated(this%dofs)) then
274 deallocate(this%dofs)
275 end if
276 if (allocated(this%interp_r)) then
277 deallocate(this%interp_r)
278 end if
279 if (allocated(this%interp_p)) then
280 deallocate(this%interp_p)
281 end if
282 end subroutine node_free
283
290 recursive subroutine tamg_matvec(this, vec_out, vec_in, lvl_out)
291 class(tamg_hierarchy_t), intent(inout) :: this
292 real(kind=rp), intent(inout) :: vec_out(:)
293 real(kind=rp), intent(inout) :: vec_in(:)
294 integer, intent(in) :: lvl_out
295 integer :: i, n, e
296 !call this%matvec_impl(vec_out, vec_in, this%nlvls, lvl_out)
297 call tamg_matvec_flat_impl(this, vec_out, vec_in, this%nlvls, lvl_out)
298 end subroutine tamg_matvec
299
307 recursive subroutine tamg_matvec_impl(this, vec_out, vec_in, lvl, lvl_out)
308 class(tamg_hierarchy_t), intent(inout) :: this
309 real(kind=rp), intent(inout) :: vec_out(:)
310 real(kind=rp), intent(inout) :: vec_in(:)
311 integer, intent(in) :: lvl
312 integer, intent(in) :: lvl_out
313 integer :: i, n, e
314
315 if (lvl .eq. 0) then
317 n = size(vec_in)
319 call this%gs_h%op(vec_in, n, gs_op_add)
320 call col2( vec_in, this%coef%mult, n)
321
322 call this%ax%compute(vec_out, vec_in, this%coef, this%msh, this%Xh)
323 call this%gs_h%op(vec_out, n, gs_op_add)
324 call this%blst%apply(vec_out, n)
325
326 if (lvl_out .ne. 0) then
327 call col2(vec_out, this%coef%mult, n)
328 end if
330 else
331 if (lvl_out .ge. lvl) then
334 associate( wrk_in => this%lvl(lvl)%wrk_in, wrk_out => this%lvl(lvl)%wrk_out)
335 n = this%lvl(lvl)%fine_lvl_dofs
336 call rzero(wrk_in, n)
337 call rzero(vec_out, this%lvl(lvl)%nnodes)
338 do n = 1, this%lvl(lvl)%nnodes
339 associate(node => this%lvl(lvl)%nodes(n))
340 do i = 1, node%ndofs
341 wrk_in( node%dofs(i) ) = wrk_in( node%dofs(i) ) + vec_in( node%gid ) * node%interp_p( i )
342 end do
343 end associate
344 end do
345
346 call this%matvec_impl(wrk_out, wrk_in, lvl-1, lvl_out)
347
349 do n = 1, this%lvl(lvl)%nnodes
350 associate(node => this%lvl(lvl)%nodes(n))
351 do i = 1, node%ndofs
352 vec_out( node%gid ) = vec_out(node%gid ) + wrk_out( node%dofs(i) ) * node%interp_r( i )
353 end do
354 end associate
355 end do
356 end associate
357 else if (lvl_out .lt. lvl) then
359 call this%matvec_impl(vec_out, vec_in, lvl-1, lvl_out)
360 else
361 call neko_error("TAMG: matvec level numbering problem.")
362 end if
363 end if
364 end subroutine tamg_matvec_impl
365
366
368 recursive subroutine tamg_matvec_flat_impl(this, vec_out, vec_in, lvl_blah, lvl_out)
369 class(tamg_hierarchy_t), intent(inout) :: this
370 real(kind=rp), intent(inout) :: vec_out(:)
371 real(kind=rp), intent(inout) :: vec_in(:)
372 integer, intent(in) :: lvl_blah
373 integer, intent(in) :: lvl_out
374 integer :: i, n, cdof, lvl
375
376 lvl = lvl_out
377 n = this%lvl(1)%fine_lvl_dofs
378 if (lvl .eq. 0) then
379 call this%ax%compute(vec_out, vec_in, this%coef, this%msh, this%Xh)
380 call this%gs_h%op(vec_out, n, gs_op_add)
381 call this%blst%apply(vec_out, n)
382 else
383 associate( wrk_in => this%lvl(1)%wrk_in, wrk_out => this%lvl(1)%wrk_out)
385 do i = 1, n
386 cdof = this%lvl(lvl)%map_finest2lvl(i)
387 wrk_in(i) = vec_in( cdof )
388 end do
389
391 call this%gs_h%op(wrk_in, n, gs_op_add)
392 call col2( wrk_in, this%coef%mult, n)
393 call this%blst%apply(wrk_in, n)
394
396 call this%ax%compute(wrk_out, wrk_in, this%coef, this%msh, this%Xh)
397 call this%gs_h%op(wrk_out, n, gs_op_add)
398 call this%blst%apply(wrk_out, n)
399
400 call col2(wrk_out, this%coef%mult, n)
401
403 call rzero(vec_out, this%lvl(lvl)%nnodes)
404 !$omp parallel do private(cdof)
405 do i = 1, n
406 cdof = this%lvl(lvl)%map_finest2lvl(i)
407 !$omp atomic
408 vec_out(cdof) = vec_out(cdof) + wrk_out( i )
409 end do
410 !$omp end parallel do
411 end associate
412 end if
413 end subroutine tamg_matvec_flat_impl
414
415
416
421 subroutine tamg_restriction_operator(this, vec_out, vec_in, lvl)
422 class(tamg_hierarchy_t), intent(inout) :: this
423 real(kind=rp), intent(inout) :: vec_out(:)
424 real(kind=rp), intent(inout) :: vec_in(:)
425 integer, intent(in) :: lvl
426 integer :: i, n, node_start, node_end, node_id
427
428 vec_out = 0d0
429 if (lvl-1 .eq. 0) then
430 call col2(vec_in, this%coef%mult, this%lvl(lvl)%fine_lvl_dofs)
431 end if
432 do n = 1, this%lvl(lvl)%nnodes
433 associate(node => this%lvl(lvl)%nodes(n))
434 do i = 1, node%ndofs
435 vec_out( node%gid ) = vec_out( node%gid ) + vec_in( node%dofs(i) ) * node%interp_r( i )
436 end do
437 end associate
438 end do
439 end subroutine tamg_restriction_operator
440
445 subroutine tamg_prolongation_operator(this, vec_out, vec_in, lvl)
446 class(tamg_hierarchy_t), intent(inout) :: this
447 real(kind=rp), intent(inout) :: vec_out(:)
448 real(kind=rp), intent(inout) :: vec_in(:)
449 integer, intent(in) :: lvl
450 integer :: i, n, node_start, node_end, node_id
451
452 vec_out = 0d0
453 do n = 1, this%lvl(lvl)%nnodes
454 associate(node => this%lvl(lvl)%nodes(n))
455 do i = 1, node%ndofs
456 vec_out( node%dofs(i) ) = vec_out( node%dofs(i) ) + vec_in( node%gid ) * node%interp_p( i )
457 end do
458 end associate
459 end do
460 if (lvl-1 .eq. 0) then
461 call this%gs_h%op(vec_out, this%lvl(lvl)%fine_lvl_dofs, gs_op_add)
462 call col2(vec_out, this%coef%mult, this%lvl(lvl)%fine_lvl_dofs)
463 call this%blst%apply(vec_out, this%lvl(lvl)%fine_lvl_dofs)
464 end if
465 end subroutine tamg_prolongation_operator
466
467
468 subroutine tamg_device_matvec_flat_impl(this, vec_out, vec_in, vec_out_d, vec_in_d, lvl_out)
469 class(tamg_hierarchy_t), intent(inout) :: this
470 real(kind=rp), intent(inout) :: vec_out(:)
471 real(kind=rp), intent(inout) :: vec_in(:)
472 type(c_ptr) :: vec_out_d
473 type(c_ptr) :: vec_in_d
474 integer, intent(in) :: lvl_out
475 integer :: i, n, cdof, lvl
476
477 lvl = lvl_out
478 n = this%lvl(1)%fine_lvl_dofs
479 if (lvl .eq. 0) then
480 call this%ax%compute(vec_out, vec_in, this%coef, this%msh, this%Xh)
481 call this%gs_h%op(vec_out, n, gs_op_add, glb_cmd_event)
482 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
483 call this%blst%apply(vec_out, n)
484 else
485
486 associate( wrk_in_d => this%lvl(1)%wrk_in_d, wrk_out_d => this%lvl(1)%wrk_out_d)
488 call device_masked_gather_copy_0(wrk_in_d, vec_in_d, this%lvl(lvl)%map_finest2lvl_d, this%lvl(lvl)%nnodes, n)
490 call this%gs_h%op(this%lvl(1)%wrk_in, n, gs_op_add, glb_cmd_event)
491 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
492 call device_col2( wrk_in_d, this%coef%mult_d, n)
493 call this%blst%apply(this%lvl(1)%wrk_in, n)
494
496 call this%ax%compute(this%lvl(1)%wrk_out, this%lvl(1)%wrk_in, this%coef, this%msh, this%Xh)
497 call this%gs_h%op(this%lvl(1)%wrk_out, n, gs_op_add, glb_cmd_event)
498 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
499 call this%blst%apply(this%lvl(1)%wrk_out, n)
500
501 call device_col2( wrk_out_d, this%coef%mult_d, n)
502
504 call device_rzero(vec_out_d, this%lvl(lvl)%nnodes)
505 call device_masked_atomic_reduction_0(vec_out_d, wrk_out_d, this%lvl(lvl)%map_finest2lvl_d, this%lvl(lvl)%nnodes, n)
506 end associate
507
508 end if
509 end subroutine tamg_device_matvec_flat_impl
510
511 subroutine tamg_device_restriction_operator(this, vec_out_d, vec_in_d, lvl)
512 class(tamg_hierarchy_t), intent(inout) :: this
513 type(c_ptr) :: vec_out_d
514 type(c_ptr) :: vec_in_d
515 integer, intent(in) :: lvl
516 integer :: i, n, m
517 n = this%lvl(lvl)%nnodes
518 m = this%lvl(lvl)%fine_lvl_dofs
519 if (lvl-1 .eq. 0) then
520 call device_col2(vec_in_d, this%coef%mult_d, m)
521 end if
522 call device_rzero(vec_out_d, n)
523 call device_masked_atomic_reduction_0(vec_out_d, vec_in_d, this%lvl(lvl)%map_f2c_d, n, m)
525
526 subroutine tamg_device_prolongation_operator(this, vec_out_d, vec_in_d, lvl, vec_out)
527 class(tamg_hierarchy_t), intent(inout) :: this
528 real(kind=rp), intent(inout) :: vec_out(:)
529 type(c_ptr) :: vec_out_d
530 type(c_ptr) :: vec_in_d
531 integer, intent(in) :: lvl
532 integer :: i, n, m
533 n = this%lvl(lvl)%nnodes
534 m = this%lvl(lvl)%fine_lvl_dofs
535 call device_masked_gather_copy_0(vec_out_d, vec_in_d, this%lvl(lvl)%map_f2c_d, n, m)
536 if (lvl-1 .eq. 0) then
537 call this%gs_h%op(vec_out, m, gs_op_add, glb_cmd_event)
538 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
539 call device_col2( vec_out_d, this%coef%mult_d, m)
540 call this%blst%apply( vec_out, m)
541 end if
543
544end module tree_amg
Map a Fortran array to a device (allocate and associate)
Definition device.F90:78
Unmap a Fortran array from a device (deassociate and free)
Definition device.F90:84
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
subroutine, public device_masked_atomic_reduction_0(a_d, b_d, mask_d, n, n_mask, strm)
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
subroutine, public device_masked_gather_copy_0(a_d, b_d, mask_d, n, n_mask, strm)
Gather a masked vector .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
subroutine, public device_stream_wait_event(stream, event, flags)
Synchronize a device stream with an event.
Definition device.F90:1471
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:52
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:63
Gather-scatter.
Definition math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:1044
subroutine, public rzero(a, n)
Zero a real vector.
Definition math.f90:233
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 the base type for TreeAMG hierarchy structure.
Definition tree_amg.f90:34
subroutine tamg_device_matvec_flat_impl(this, vec_out, vec_in, vec_out_d, vec_in_d, lvl_out)
Definition tree_amg.f90:469
subroutine tamg_device_prolongation_operator(this, vec_out_d, vec_in_d, lvl, vec_out)
Definition tree_amg.f90:527
subroutine lvl_free(this)
deallocate tamg level
Definition tree_amg.f90:212
recursive subroutine tamg_matvec_impl(this, vec_out, vec_in, lvl, lvl_out)
Matrix vector product using the TreeAMG hierarchy structure b=Ax done as vec_out = A * vec_in This is...
Definition tree_amg.f90:308
subroutine, public tamg_node_init(node, gid, ndofs)
Initialization of a TreeAMG tree node.
Definition tree_amg.f90:255
subroutine tamg_restriction_operator(this, vec_out, vec_in, lvl)
Restriction operator for TreeAMG. vec_out = R * vec_in.
Definition tree_amg.f90:422
subroutine tamg_prolongation_operator(this, vec_out, vec_in, lvl)
Prolongation operator for TreeAMG. vec_out = P * vec_in.
Definition tree_amg.f90:446
subroutine node_free(this)
deallocate tamg tree node
Definition tree_amg.f90:271
subroutine, public tamg_lvl_init(tamg_lvl, lvl, nnodes, ndofs)
Initialization of a TreeAMG level.
Definition tree_amg.f90:186
subroutine tamg_free(this)
deallocate tamg hierarchy
Definition tree_amg.f90:161
recursive subroutine tamg_matvec(this, vec_out, vec_in, lvl_out)
Wrapper for matrix vector product using the TreeAMG hierarchy structure b=Ax done as vec_out = A * ve...
Definition tree_amg.f90:291
subroutine tamg_device_restriction_operator(this, vec_out_d, vec_in_d, lvl)
Definition tree_amg.f90:512
subroutine tamg_init(this, ax, xh, coef, msh, gs_h, nlvls, blst)
Initialization of TreeAMG hierarchy.
Definition tree_amg.f90:126
recursive subroutine tamg_matvec_flat_impl(this, vec_out, vec_in, lvl_blah, lvl_out)
Ignore this. For piecewise constant, can create index map directly to finest level.
Definition tree_amg.f90:369
Utilities.
Definition utils.f90:35
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:49
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
The function space for the SEM solution fields.
Definition space.f90:63
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:87
Type for storing TreeAMG level information.
Definition tree_amg.f90:68
Type for storing TreeAMG tree node information.
Definition tree_amg.f90:54