Neko 1.99.2
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, 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
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 (neko_bcknd_device .eq. 1) then
216 call device_deassociate(this%wrk_in)
217 call device_deassociate(this%wrk_out)
218 call device_deassociate(this%map_f2c)
219 call device_deassociate(this%map_finest2lvl)
220 end if
221 if (allocated(this%nodes)) then
222 do i = 1, this%nnodes
223 call this%nodes(i)%free()
224 end do
225 deallocate(this%nodes)
226 end if
227 if (allocated(this%wrk_in)) then
228 deallocate(this%wrk_in)
229 end if
230 if (allocated(this%wrk_out)) then
231 deallocate(this%wrk_out)
232 end if
233 if (allocated(this%map_f2c)) then
234 deallocate(this%map_f2c)
235 end if
236 if (allocated(this%map_finest2lvl)) then
237 deallocate(this%map_finest2lvl)
238 end if
239 if (neko_bcknd_device .eq. 1) then
240 call device_free(this%wrk_in_d)
241 call device_free(this%wrk_out_d)
242 call device_free(this%map_f2c_d)
243 call device_free(this%map_finest2lvl_d)
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 do i = 1, n
405 cdof = this%lvl(lvl)%map_finest2lvl(i)
406 vec_out(cdof) = vec_out(cdof) + wrk_out( i )
407 end do
408 end associate
409 end if
410 end subroutine tamg_matvec_flat_impl
411
412
413
418 subroutine tamg_restriction_operator(this, vec_out, vec_in, lvl)
419 class(tamg_hierarchy_t), intent(inout) :: this
420 real(kind=rp), intent(inout) :: vec_out(:)
421 real(kind=rp), intent(inout) :: vec_in(:)
422 integer, intent(in) :: lvl
423 integer :: i, n, node_start, node_end, node_id
424
425 vec_out = 0d0
426 if (lvl-1 .eq. 0) then
427 call col2(vec_in, this%coef%mult, this%lvl(lvl)%fine_lvl_dofs)
428 end if
429 do n = 1, this%lvl(lvl)%nnodes
430 associate(node => this%lvl(lvl)%nodes(n))
431 do i = 1, node%ndofs
432 vec_out( node%gid ) = vec_out( node%gid ) + vec_in( node%dofs(i) ) * node%interp_r( i )
433 end do
434 end associate
435 end do
436 end subroutine tamg_restriction_operator
437
442 subroutine tamg_prolongation_operator(this, vec_out, vec_in, lvl)
443 class(tamg_hierarchy_t), intent(inout) :: this
444 real(kind=rp), intent(inout) :: vec_out(:)
445 real(kind=rp), intent(inout) :: vec_in(:)
446 integer, intent(in) :: lvl
447 integer :: i, n, node_start, node_end, node_id
448
449 vec_out = 0d0
450 do n = 1, this%lvl(lvl)%nnodes
451 associate(node => this%lvl(lvl)%nodes(n))
452 do i = 1, node%ndofs
453 vec_out( node%dofs(i) ) = vec_out( node%dofs(i) ) + vec_in( node%gid ) * node%interp_p( i )
454 end do
455 end associate
456 end do
457 if (lvl-1 .eq. 0) then
458 call this%gs_h%op(vec_out, this%lvl(lvl)%fine_lvl_dofs, gs_op_add)
459 call col2(vec_out, this%coef%mult, this%lvl(lvl)%fine_lvl_dofs)
460 call this%blst%apply(vec_out, n)
461 end if
462 end subroutine tamg_prolongation_operator
463
464
465 subroutine tamg_device_matvec_flat_impl(this, vec_out, vec_in, vec_out_d, vec_in_d, lvl_out)
466 class(tamg_hierarchy_t), intent(inout) :: this
467 real(kind=rp), intent(inout) :: vec_out(:)
468 real(kind=rp), intent(inout) :: vec_in(:)
469 type(c_ptr) :: vec_out_d
470 type(c_ptr) :: vec_in_d
471 integer, intent(in) :: lvl_out
472 integer :: i, n, cdof, lvl
473
474 lvl = lvl_out
475 n = this%lvl(1)%fine_lvl_dofs
476 if (lvl .eq. 0) then
477 call this%ax%compute(vec_out, vec_in, this%coef, this%msh, this%Xh)
478 call this%gs_h%op(vec_out, n, gs_op_add, glb_cmd_event)
479 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
480 call this%blst%apply(vec_out, n)
481 else
482
483 associate( wrk_in_d => this%lvl(1)%wrk_in_d, wrk_out_d => this%lvl(1)%wrk_out_d)
485 call device_masked_gather_copy_0(wrk_in_d, vec_in_d, this%lvl(lvl)%map_finest2lvl_d, this%lvl(lvl)%nnodes, n)
487 call this%gs_h%op(this%lvl(1)%wrk_in, n, gs_op_add, glb_cmd_event)
488 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
489 call device_col2( wrk_in_d, this%coef%mult_d, n)
490 call this%blst%apply(this%lvl(1)%wrk_in, n)
491
493 call this%ax%compute(this%lvl(1)%wrk_out, this%lvl(1)%wrk_in, this%coef, this%msh, this%Xh)
494 call this%gs_h%op(this%lvl(1)%wrk_out, n, gs_op_add, glb_cmd_event)
495 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
496 call this%blst%apply(this%lvl(1)%wrk_out, n)
497
498 call device_col2( wrk_out_d, this%coef%mult_d, n)
499
501 call device_rzero(vec_out_d, this%lvl(lvl)%nnodes)
502 call device_masked_atomic_reduction_0(vec_out_d, wrk_out_d, this%lvl(lvl)%map_finest2lvl_d, this%lvl(lvl)%nnodes, n)
503 end associate
504
505 end if
506 end subroutine tamg_device_matvec_flat_impl
507
508 subroutine tamg_device_restriction_operator(this, vec_out_d, vec_in_d, lvl)
509 class(tamg_hierarchy_t), intent(inout) :: this
510 type(c_ptr) :: vec_out_d
511 type(c_ptr) :: vec_in_d
512 integer, intent(in) :: lvl
513 integer :: i, n, m
514 n = this%lvl(lvl)%nnodes
515 m = this%lvl(lvl)%fine_lvl_dofs
516 if (lvl-1 .eq. 0) then
517 call device_col2(vec_in_d, this%coef%mult_d, m)
518 end if
519 call device_rzero(vec_out_d, n)
520 call device_masked_atomic_reduction_0(vec_out_d, vec_in_d, this%lvl(lvl)%map_f2c_d, n, m)
522
523 subroutine tamg_device_prolongation_operator(this, vec_out_d, vec_in_d, lvl, vec_out)
524 class(tamg_hierarchy_t), intent(inout) :: this
525 real(kind=rp), intent(inout) :: vec_out(:)
526 type(c_ptr) :: vec_out_d
527 type(c_ptr) :: vec_in_d
528 integer, intent(in) :: lvl
529 integer :: i, n, m
530 n = this%lvl(lvl)%nnodes
531 m = this%lvl(lvl)%fine_lvl_dofs
532 call device_masked_gather_copy_0(vec_out_d, vec_in_d, this%lvl(lvl)%map_f2c_d, n, m)
533 if (lvl-1 .eq. 0) then
534 call this%gs_h%op(vec_out, m, gs_op_add, glb_cmd_event)
535 call device_stream_wait_event(glb_cmd_queue, glb_cmd_event, 0)
536 call device_col2( vec_out_d, this%coef%mult_d, m)
537 call this%blst%apply( vec_out, n)
538 end if
540
541end module tree_amg
Deassociate a Fortran array from a device pointer.
Definition device.F90:95
Map a Fortran array to a device (allocate and associate)
Definition device.F90:77
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_free(x_d)
Deallocate memory on the device.
Definition device.F90:219
subroutine, public device_stream_wait_event(stream, event, flags)
Synchronize a device stream with an event.
Definition device.F90:1208
type(c_ptr), bind(C), public glb_cmd_queue
Global command queue.
Definition device.F90:51
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:62
Gather-scatter.
Definition math.f90:60
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 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:466
subroutine tamg_device_prolongation_operator(this, vec_out_d, vec_in_d, lvl, vec_out)
Definition tree_amg.f90:524
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:419
subroutine tamg_prolongation_operator(this, vec_out, vec_in, lvl)
Prolongation operator for TreeAMG. vec_out = P * vec_in.
Definition tree_amg.f90:443
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:509
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:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
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