Loading [MathJax]/extensions/tex2jax.js
Neko 0.9.99
A portable framework for high-order spectral element flow simulations
All Classes Namespaces Files Functions Variables Typedefs Enumerator Macros Pages
tree_amg_aggregate.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!
36 use utils, only : neko_error, linear_index
37 use num_types, only : rp, dp, i8
38 use comm
39 use mesh, only : mesh_t
40 use logger, only : neko_log, log_size
41 implicit none
42
43 type, public :: tamg_agg_monitor_t
45 integer :: level
46 integer :: num_dofs
47 integer :: num_aggs
49 integer :: phase1_naggs
50 integer :: phase1_ndof_aggregated
51 integer :: phase2_naggs
52 integer :: phase2_ndof_aggregated
53 end type tamg_agg_monitor_t
54
55contains
56
64 subroutine aggregate_finest_level(tamg, lx, ly, lz, ne)
65 type(tamg_hierarchy_t), intent(inout) :: tamg
66 integer, intent(in) :: lx, ly, lz, ne
67 integer :: i, j, k, l, nl, nt
68 integer :: lid, gid_ptr
69 integer :: lvl_id
70 lvl_id = 1
72 nl = lx*ly*lz
73 nt = nl*ne
75 call tamg_lvl_init(tamg%lvl(lvl_id), lvl_id, ne, nt)
76 gid_ptr = 1
77 do l = 1, ne
78 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
79 tamg%lvl(lvl_id)%nodes_gid(l) = l
80 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, nl)
82 lid = 0
83 do k = 1, lz
84 do j = 1, ly
85 do i = 1, lx
86 lid = lid + 1
87 tamg%lvl(lvl_id)%nodes(l)%dofs(lid) = linear_index(i,j,k,l,lx,ly,lz)
88
89 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = linear_index(i,j,k,l,lx,ly,lz)
90 !tamg%lvl(lvl_id)%nodes_gids(gid_ptr) = l
91 tamg%lvl(lvl_id)%map_f2c(linear_index(i,j,k,l,lx,ly,lz)) = l
92 gid_ptr = gid_ptr + 1
93 end do
94 end do
95 end do
96 end do
97 tamg%lvl(lvl_id)%nodes_ptr(ne+1) = gid_ptr
98
99 call aggregation_monitor_finest(lvl_id,nt,ne)
100
101 end subroutine aggregate_finest_level
102
112 subroutine agg_greedy_first_pass(naggs, max_aggs, n_elements, &
113 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
114 integer, intent(inout):: naggs
115 integer, intent(in) :: max_aggs, n_elements
116 integer, intent(in) :: facet_neigh(:,:)
117 integer, intent(in) :: offset_el, n_facet
118 integer, intent(inout) :: is_aggregated(:)
119 integer, allocatable, intent(inout) :: aggregate_size(:)
120 integer, allocatable :: as_tmp(:)
121 real(kind=dp) :: random_value
122 integer :: i, side, nhbr
123 logical :: no_nhbr_agg
124
125!-----! do while (naggs .le. max_aggs)
126!-----! call random_number(random_value)
127!-----! i = floor(random_value * n_elements + 1)
128 do i = 1, n_elements
129 if (is_aggregated(i) .eq. -1) then
131 no_nhbr_agg = .true.
132 do side = 1, n_facet
133 nhbr = facet_neigh(side, i) - offset_el
134 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
135 if (is_aggregated(nhbr) .ne. -1) then
136 no_nhbr_agg = .false.
137 end if
138 end if
139 end do! side
140
142 if (no_nhbr_agg) then
143 naggs = naggs + 1
144 is_aggregated(i) = naggs
145 if(size(aggregate_size).lt.naggs) then
146 allocate(as_tmp(naggs + 20))
147 as_tmp(1:size(aggregate_size)) = aggregate_size
148 call move_alloc(as_tmp, aggregate_size)
149 end if
150 aggregate_size(naggs) = 1
151 do side = 1, n_facet
152 nhbr = facet_neigh(side, i) - offset_el
153 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
154 if (is_aggregated(nhbr) .eq. -1) then
155 is_aggregated(nhbr) = naggs
156 aggregate_size(naggs) = aggregate_size(naggs) + 1
157 end if
158 end if
159 end do! side
160 end if! no_nhbr_agg
161 end if! is_aggregated(i)
162 end do
163 end subroutine agg_greedy_first_pass
164
175 subroutine agg_greedy_second_pass(naggs, max_aggs, n_elements, &
176 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
177 integer, intent(inout):: naggs
178 integer, intent(in) :: max_aggs, n_elements
179 integer, intent(in) :: facet_neigh(:,:)
180 integer, intent(in) :: offset_el, n_facet
181 integer, intent(inout) :: is_aggregated(:)
182 integer, intent(inout) :: aggregate_size(:)
183 integer :: i, side, nhbr
184 integer :: tnt_agg, tst_agg, tnt_size, tst_size
185
187 do i = 1, n_elements
188 if (is_aggregated(i) .eq. -1) then
190 tnt_agg = -1
191 tnt_size = 999!TODO: replace with large number
192 tst_agg = -1
193 tst_size = 999!TODO: replace with large number
194 do side = 1, n_facet
195 nhbr = facet_neigh(side, i) - offset_el
196 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
197 if (is_aggregated(nhbr) .ne. -1) then
198 tst_agg = is_aggregated(nhbr)
199 tst_size = aggregate_size(tst_agg)
200 if (tst_size .lt. tnt_size) then
201 tnt_size = tst_size
202 tnt_agg = tst_agg
203 end if
204 end if
205 end if
206 end do
207
208 if (tnt_agg .ne. -1) then
210 is_aggregated(i) = tnt_agg
211 aggregate_size(tnt_agg) = aggregate_size(tnt_agg) + 1
212 else
214 naggs = naggs + 1
215 if (naggs .gt. size(aggregate_size)) then!TODO: another movealoc here? Error? the max_aggs needs to change though...
216 call neko_error("Creating too many aggregates... something might be wrong... try increasing max_aggs")
217 end if
218 is_aggregated(i) = naggs
219 aggregate_size(naggs) = 1
221 do side = 1, n_facet
222 nhbr = facet_neigh(side, i) - offset_el
223 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
224 if (is_aggregated(nhbr) .eq. -1) then
225 aggregate_size(naggs) = aggregate_size(naggs) + 1
226 is_aggregated(nhbr) = naggs
227 end if
228 end if
229 end do
230 end if
231
232 end if
233 end do
234 end subroutine agg_greedy_second_pass
235
245 subroutine agg_fill_nhbr_info(agg_nhbr, n_agg_nhbr, n_elements, &
246 facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
247 integer, intent(inout) :: agg_nhbr(:,:)
248 integer, intent(inout) :: n_agg_nhbr
249 integer, intent(in) :: n_elements
250 integer, intent(in) :: facet_neigh(:,:)
251 integer, intent(in) :: offset_el, n_facet
252 integer, intent(inout) :: is_aggregated(:)
253 integer, intent(inout) :: aggregate_size(:)
254 integer :: i, j, side, nhbr, tst_agg, tnt_agg
255 logical :: agg_added
256
257 n_agg_nhbr = 0
258
259 do i = 1, n_elements!TODO: this is the lazy expensive way...
260 tnt_agg = is_aggregated(i)
261 do side = 1, n_facet
262 nhbr = facet_neigh(side,i) - offset_el
263 if ((nhbr .gt. 0).and.(nhbr .le. n_elements)) then
264 tst_agg = is_aggregated(nhbr)
265 if (tst_agg .le. 0) call neko_error("Unaggregated element detected. We do not want to handle that here...")
266 if (tst_agg .ne. tnt_agg) then
267 agg_added = .false.
268 do j = 1, 50!TODO: this hard-coded value
269 if ((agg_nhbr(j,tnt_agg) .eq. tst_agg)) then
270 agg_added = .true.
271 else if ((agg_nhbr(j,tnt_agg).eq.-1).and.(.not.agg_added)) then
272 agg_nhbr(j,tnt_agg) = tst_agg
273 agg_added = .true.
274 n_agg_nhbr = max(n_agg_nhbr, j)
275 end if
276 end do! j
277 if (.not.agg_added) call neko_error("Aggregates have too many neighbors... probably. Or some other error.")
278 end if
279 end if
280 end do! side
281 end do! i
282 end subroutine agg_fill_nhbr_info
283
290 subroutine aggregate_greedy(tamg, lvl_id, max_aggs, facet_neigh, agg_nhbr)
291 type(tamg_hierarchy_t), intent(inout) :: tamg
292 integer, intent(in) :: lvl_id
293 integer, intent(in) :: max_aggs
294 integer, intent(in) :: facet_neigh(:,:)
295 integer, intent(inout), allocatable :: agg_nhbr(:,:)
296 integer, allocatable :: is_aggregated(:)
297 integer, allocatable :: aggregate_size(:)
298 integer :: n_elements, naggs, n_facet, offset_el
299 integer :: i, j, l, ntot, n_agg_facet, gid_ptr
300
301 if (lvl_id .lt. 2) then
302 call neko_error("For now, can only use greedy agg after elms have been aggregated to points (level 1)")
303 else if (lvl_id .eq. 2) then
304 n_facet = 6
305 offset_el = tamg%msh%offset_el
306 else
307 n_facet = 50!TODO: this hard-coded value. how many neighbors can there be?
308 offset_el = 0
309 end if
310
311 naggs = 0
312 n_elements = tamg%lvl(lvl_id-1)%nnodes
313 allocate( is_aggregated( n_elements ) )
314 allocate( aggregate_size( max_aggs ) )
315
317 is_aggregated = -1
318
320 aggregate_size = huge(i)!999999
321
323 call agg_greedy_first_pass(naggs, max_aggs, n_elements, &
324 facet_neigh, offset_el, n_facet, &
325 is_aggregated, aggregate_size)
326
327 call aggregation_monitor_phase1(lvl_id, n_elements, naggs, is_aggregated)
328
330 call agg_greedy_second_pass(naggs, max_aggs, n_elements, &
331 facet_neigh, offset_el, n_facet, &
332 is_aggregated, aggregate_size)
333
334 call aggregation_monitor_phase2(lvl_id, n_elements, naggs, is_aggregated)
335
336 if (.true.) then
337 allocate( agg_nhbr(50, max_aggs*2) )!TODO: this hard-coded n_facet value (20)...
338 agg_nhbr = -1
339 call agg_fill_nhbr_info(agg_nhbr, n_agg_facet, n_elements, &
340 facet_neigh, offset_el, n_facet, &
341 is_aggregated, aggregate_size)
342 end if
343
345 ntot = 0
346 do l = 1, naggs
347 ntot = ntot + aggregate_size(l)
348 end do
350 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, naggs, ntot)
351 ntot = 0
352 gid_ptr = 1
353 do l = 1, naggs
354 tamg%lvl(lvl_id)%nodes_ptr(l) = gid_ptr
355 tamg%lvl(lvl_id)%nodes_gid(l) = l
356 call tamg_node_init( tamg%lvl(lvl_id)%nodes(l), l, aggregate_size(l))
357 j = 0
358 do i = 1, n_elements!TODO: this is the lazy expensive way...
359 if (is_aggregated(i) .eq. l) then
360 j = j+1
361 tamg%lvl(lvl_id)%nodes(l)%dofs(j) = i
362
363 tamg%lvl(lvl_id)%nodes_dofs(gid_ptr) = i
364 !tamg%lvl(lvl_id)%nodes_gids(gid_ptr) = l
365 tamg%lvl(lvl_id)%map_f2c(i) = l
366 gid_ptr = gid_ptr + 1
367 end if
368 end do
369 if (j .ne. tamg%lvl(lvl_id)%nodes(l)%ndofs) then
370 call neko_error("Aggregation problem. Not enough dofs in node.")
371 end if
372 ntot = ntot + aggregate_size(l)
373 end do
374 tamg%lvl(lvl_id)%nodes_ptr(naggs+1) = gid_ptr
375
376 call aggregation_monitor_final(lvl_id,ntot,naggs)
377
378 deallocate( is_aggregated )
379 deallocate( aggregate_size )
380 end subroutine aggregate_greedy
381
385 subroutine aggregate_end( tamg, lvl_id)
386 type(tamg_hierarchy_t), intent(inout) :: tamg
387 integer, intent(in) :: lvl_id
388 integer :: nt, i
390 nt = tamg%lvl(lvl_id-1)%nnodes
392 call tamg_lvl_init( tamg%lvl(lvl_id), lvl_id, 1, nt)
393
395 call tamg_node_init( tamg%lvl(lvl_id)%nodes(1), 1, nt)
397 do i = 1, tamg%lvl(lvl_id-1)%nnodes
398 tamg%lvl(lvl_id)%nodes(1)%dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
399
400 tamg%lvl(lvl_id)%nodes_dofs(i) = tamg%lvl(lvl_id-1)%nodes(i)%gid
401 tamg%lvl(lvl_id)%map_f2c(i) = 1
402 end do
403
404 tamg%lvl(lvl_id)%nodes_ptr(1) = 1
405 tamg%lvl(lvl_id)%nodes_ptr(2) = 2
406 tamg%lvl(lvl_id)%nodes_gid(1) = 1
407 end subroutine aggregate_end
408
409 subroutine aggregation_monitor_finest(lvl,ndof,nagg)
410 integer, intent(in) :: lvl,ndof,nagg
411 integer :: na_max, na_min, na_avg, na_sum
412 character(len=LOG_SIZE) :: log_buf
413
414 write(log_buf, '(A8,I2,A37)') '-- level',lvl,'-- Aggregation: Element-as-Aggregate'
415 !write(log_buf, '(A44)') 'Aggregation: Element-as-Aggregate'
416 call neko_log%message(log_buf)
417
418 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
419 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
420 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
421 na_avg = na_sum / pe_size
422 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
423 call neko_log%message(log_buf)
424
425 end subroutine aggregation_monitor_finest
426
427 subroutine aggregation_monitor_phase1(lvl,ndof,nagg,is_aggregated)
428 integer, intent(in) :: lvl,ndof,nagg
429 integer, intent(in) :: is_aggregated(:)
430 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
431 real(kind=rp) :: agg_frac
432 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
433 character(len=LOG_SIZE) :: log_buf
434 num_aggregated = 0
435 do i = 1, ndof
436 if (is_aggregated(i) .gt. -1) then
437 num_aggregated = num_aggregated + 1
438 end if
439 end do
440
441 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase1'
442 write(log_buf, '(A27)') 'Aggregation: phase1'
443 call neko_log%message(log_buf)
444
445 loc_aggd = int(num_aggregated, i8)
446 loc_dof = int(ndof, i8)
447 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
448 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
449 agg_frac = real(glb_aggd,rp) / real(glb_dof,rp)
450
451 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
452 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
453 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
454 na_avg = na_sum / pe_size
455 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
456 call neko_log%message(log_buf)
457
458 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
459 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
460 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
461 na_avg = na_sum / pe_size
462 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (',na_min,',',na_avg,',',na_max,')',(agg_frac*100)
463 call neko_log%message(log_buf)
464 end subroutine aggregation_monitor_phase1
465
466 subroutine aggregation_monitor_phase2(lvl,ndof,nagg,is_aggregated)
467 integer, intent(in) :: lvl,ndof,nagg
468 integer, intent(in) :: is_aggregated(:)
469 integer :: num_aggregated, i, na_max, na_min, na_avg, na_sum
470 real(kind=rp) :: agg_frac
471 integer(kind=i8) :: glb_dof, glb_aggd, loc_dof, loc_aggd
472 character(len=LOG_SIZE) :: log_buf
473 num_aggregated = 0
474 do i = 1, ndof
475 if (is_aggregated(i) .gt. -1) then
476 num_aggregated = num_aggregated + 1
477 end if
478 end do
479 !write(log_buf, '(A8,I2,A24)') '-- level',lvl,'-- Aggregation: phase2'
480 write(log_buf, '(A27)') 'Aggregation: phase2'
481 call neko_log%message(log_buf)
482
483 loc_aggd = int(num_aggregated, i8)
484 loc_dof = int(ndof, i8)
485 call mpi_allreduce(loc_dof, glb_dof, 1, mpi_integer8, mpi_sum, neko_comm)
486 call mpi_allreduce(loc_aggd, glb_aggd, 1, mpi_integer8, mpi_sum, neko_comm)
487 agg_frac = real(glb_aggd,rp) / real(glb_dof,rp)
488
489 call mpi_allreduce(nagg, na_max, 1, mpi_integer, mpi_max, neko_comm)
490 call mpi_allreduce(nagg, na_min, 1, mpi_integer, mpi_min, neko_comm)
491 call mpi_allreduce(nagg, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
492 na_avg = na_sum / pe_size
493 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1)') 'Number of Aggregates: (',na_min,',',na_avg,',',na_max,')'
494 call neko_log%message(log_buf)
495
496 call mpi_allreduce(num_aggregated, na_max, 1, mpi_integer, mpi_max, neko_comm)
497 call mpi_allreduce(num_aggregated, na_min, 1, mpi_integer, mpi_min, neko_comm)
498 call mpi_allreduce(num_aggregated, na_sum, 1, mpi_integer, mpi_sum, neko_comm)
499 na_avg = na_sum / pe_size
500 write(log_buf, '(A35,I6,A1,I6,A1,I6,A1,F6.2)') 'Aggregated: (',na_min,',',na_avg,',',na_max,')',(agg_frac*100)
501 call neko_log%message(log_buf)
502 end subroutine aggregation_monitor_phase2
503
504 subroutine aggregation_monitor_final(lvl,ndof,nagg)
505 integer, intent(in) :: lvl,ndof,nagg
506 character(len=LOG_SIZE) :: log_buf
507 !TODO: calculate min and max agg size
508 !write(log_buf, '(A8,I2,A23,I6)') '-- level',lvl,'-- Aggregation: Done.', nagg
509 write(log_buf, '(A26,I6)') 'Aggregation: Done.', nagg
510 call neko_log%message(log_buf)
511 end subroutine aggregation_monitor_final
512
513end module tree_amg_aggregate
double real
Definition comm.F90:1
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
Defines a mesh.
Definition mesh.f90:34
integer, parameter, public i8
Definition num_types.f90:7
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements an aggregation for TreeAMG hierarchy structure.
subroutine agg_greedy_first_pass(naggs, max_aggs, n_elements, facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
First pass of a greedy aggregation Loop through all dofs and aggregate on dof that has all unaggregat...
subroutine aggregate_end(tamg, lvl_id)
Aggregate all dofs to a single point to form a tree-like structure.
subroutine aggregation_monitor_phase2(lvl, ndof, nagg, is_aggregated)
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.
subroutine aggregation_monitor_final(lvl, ndof, nagg)
subroutine agg_greedy_second_pass(naggs, max_aggs, n_elements, facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
Second pass of a greedy aggregation Loop through all unaggregated dofs and add them to a neighboring ...
subroutine aggregation_monitor_finest(lvl, ndof, nagg)
subroutine agg_fill_nhbr_info(agg_nhbr, n_agg_nhbr, n_elements, facet_neigh, offset_el, n_facet, is_aggregated, aggregate_size)
Create information on which aggregates are "adjacent" to eachother.
subroutine aggregation_monitor_phase1(lvl, ndof, nagg, is_aggregated)
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
pure integer function, public linear_index(i, j, k, l, lx, ly, lz)
Compute the address of a (i,j,k,l) array with sizes (1:lx, 1:ly, 1:lz, :)
Definition utils.f90:174
Type for a TreeAMG hierarchy.
Definition tree_amg.f90:94
#define max(a, b)
Definition tensor.cu:40