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
mesh.f90
Go to the documentation of this file.
1! Copyright (c) 2018-2023, 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!
34module mesh
35 use num_types, only : rp, dp, i8
36 use point, only : point_t
37 use element, only : element_t
38 use hex
39 use quad
40 use utils, only : neko_error, neko_warning
42 use tuple, only : tuple_i4_t, tuple4_i4_t
43 use htable
44 use datadist
45 use distdata
46 use comm
48 use math
49 use uset, only : uset_i8_t
50 use curve, only : curve_t
51 use logger, only : log_size
52 implicit none
53 private
54
56 integer, public, parameter :: neko_msh_max_zlbls = 20
58 integer, public, parameter :: neko_msh_max_zlbl_len = 40
59
60 type, private :: mesh_element_t
61 class(element_t), allocatable :: e
62 end type mesh_element_t
63
64 type, public :: mesh_t
65 integer :: nelv
66 integer :: npts
67 integer :: gdim
68 integer :: mpts
69 integer :: mfcs
70 integer :: meds
71
72 integer :: glb_nelv
73 integer :: glb_mpts
74 integer :: glb_mfcs
75 integer :: glb_meds
76
77 integer :: offset_el
78 integer :: max_pts_id
79
80 type(point_t), allocatable :: points(:)
81 type(mesh_element_t), allocatable :: elements(:)
82 logical, allocatable :: dfrmd_el(:)
83
84 type(htable_i4_t) :: htp
85 type(htable_i4t4_t) :: htf
86 type(htable_i4t2_t) :: hte
87 type(htable_i4_t) :: htel
88
89
90 integer, allocatable :: facet_neigh(:,:)
91
95 class(htable_t), allocatable :: facet_map
96 type(stack_i4_t), allocatable :: point_neigh(:)
97
98 type(distdata_t) :: ddata
99 logical, allocatable :: neigh(:)
100 integer, allocatable :: neigh_order(:)
101
102 integer(2), allocatable :: facet_type(:,:)
103
104 type(facet_zone_t), allocatable :: labeled_zones(:)
105 type(facet_zone_periodic_t) :: periodic
106 type(curve_t) :: curve
107
108 logical :: lconn = .false.
109 logical :: ldist = .false.
110 logical :: lnumr = .false.
111 logical :: lgenc = .true.
112
115 procedure(mesh_deform), pass(msh), pointer :: apply_deform => null()
116 contains
117 procedure, private, pass(this) :: init_nelv => mesh_init_nelv
118 procedure, private, pass(this) :: init_dist => mesh_init_dist
119 procedure, private, pass(this) :: add_quad => mesh_add_quad
120 procedure, private, pass(this) :: add_hex => mesh_add_hex
121 procedure, private, pass(this) :: add_edge => mesh_add_edge
122 procedure, private, pass(this) :: add_face => mesh_add_face
123 procedure, private, pass(this) :: add_point => mesh_add_point
124 procedure, private, pass(this) :: get_local_point => mesh_get_local_point
125 procedure, private, pass(this) :: get_local_edge => mesh_get_local_edge
126 procedure, private, pass(this) :: get_local_facet => mesh_get_local_facet
127 procedure, private, pass(this) :: get_global_edge => mesh_get_global_edge
128 procedure, private, pass(this) :: get_global_facet => mesh_get_global_facet
129 procedure, private, pass(this) :: is_shared_point => mesh_is_shared_point
130 procedure, private, pass(this) :: is_shared_edge => mesh_is_shared_edge
131 procedure, private, pass(this) :: is_shared_facet => mesh_is_shared_facet
132 procedure, pass(this) :: free => mesh_free
133 procedure, pass(this) :: finalize => mesh_finalize
134 procedure, pass(this) :: mark_periodic_facet => mesh_mark_periodic_facet
135 procedure, pass(this) :: mark_labeled_facet => mesh_mark_labeled_facet
136 procedure, pass(this) :: mark_curve_element => mesh_mark_curve_element
137 procedure, pass(this) :: apply_periodic_facet => mesh_apply_periodic_facet
138 procedure, pass(this) :: all_deformed => mesh_all_deformed
139 procedure, pass(this) :: get_facet_ids => mesh_get_facet_ids
140 procedure, pass(this) :: reset_periodic_ids => mesh_reset_periodic_ids
141 procedure, pass(this) :: create_periodic_ids => mesh_create_periodic_ids
142 procedure, pass(this) :: generate_conn => mesh_generate_conn
143 procedure, pass(this) :: have_point_glb_idx => mesh_have_point_glb_idx
145 generic :: init => init_nelv, init_dist
147 generic :: add_element => add_quad, add_hex
150 generic :: get_local => get_local_point, get_local_edge, get_local_facet
153 generic :: get_global => get_global_edge, get_global_facet
155 generic :: is_shared => is_shared_point, is_shared_edge, is_shared_facet
156 end type mesh_t
157
158 abstract interface
159 subroutine mesh_deform(msh, x, y, z, lx, ly, lz)
160 import mesh_t
161 import rp
162 class(mesh_t) :: msh
163 integer, intent(in) :: lx, ly, lz
164 real(kind=rp), intent(inout) :: x(lx, ly, lz, msh%nelv)
165 real(kind=rp), intent(inout) :: y(lx, ly, lz, msh%nelv)
166 real(kind=rp), intent(inout) :: z(lx, ly, lz, msh%nelv)
167 end subroutine mesh_deform
168 end interface
169
170 public :: mesh_deform
171
172contains
173
175 subroutine mesh_init_nelv(this, gdim, nelv)
176 class(mesh_t), intent(inout) :: this
177 integer, intent(in) :: gdim
178 integer, intent(in) :: nelv
179 integer :: ierr
180 character(len=LOG_SIZE) :: log_buf
181
182 call this%free()
183
184 this%nelv = nelv
185 this%gdim = gdim
186
187 if (this%nelv < 1) then
188 write(log_buf, '(A,I0,A)') 'MPI rank ', pe_rank, ' has zero elements'
189 call neko_warning(log_buf)
190 end if
191
192 call mpi_allreduce(this%nelv, this%glb_nelv, 1, &
193 mpi_integer, mpi_sum, neko_comm, ierr)
194
195 this%offset_el = 0
196 call mpi_exscan(this%nelv, this%offset_el, 1, &
197 mpi_integer, mpi_sum, neko_comm, ierr)
198
199 call mesh_init_common(this)
200
201 end subroutine mesh_init_nelv
202
204 subroutine mesh_init_dist(this, gdim, dist)
205 class(mesh_t), intent(inout) :: this
206 integer, intent(in) :: gdim
207 type(linear_dist_t), intent(in) :: dist
208 character(len=LOG_SIZE) :: log_buf
209
210 call this%free()
211
212 this%nelv = dist%num_local()
213 if (this%nelv < 1) then
214 write(log_buf, '(A,I0,A)') 'MPI rank ', pe_rank, ' has zero elements'
215 call neko_warning(log_buf)
216 end if
217 this%glb_nelv = dist%num_global()
218 this%offset_el = dist%start_idx()
219 this%gdim = gdim
220
221 call mesh_init_common(this)
222
223 end subroutine mesh_init_dist
224
225 subroutine mesh_init_common(this)
226 type(mesh_t), intent(inout) :: this
227 integer :: i
228 type(tuple_i4_t) :: facet_data
229
230 this%max_pts_id = 0
231
232 allocate(this%elements(this%nelv))
233 allocate(this%dfrmd_el(this%nelv))
234 if (this%gdim .eq. 3) then
235 do i = 1, this%nelv
236 allocate(hex_t::this%elements(i)%e)
237 end do
238 this%npts = neko_hex_npts
240 if (this%lgenc) then
241 allocate(htable_i4t4_t::this%facet_map)
242 select type (fmp => this%facet_map)
243 type is(htable_i4t4_t)
244 call fmp%init(this%nelv, facet_data)
245 end select
246
247 allocate(this%facet_neigh(neko_hex_nfcs, this%nelv))
248
249 call this%htf%init(this%nelv * neko_hex_nfcs, i)
250 call this%hte%init(this%nelv * neko_hex_neds, i)
251 end if
252 else if (this%gdim .eq. 2) then
253 do i = 1, this%nelv
254 allocate(quad_t::this%elements(i)%e)
255 end do
256 this%npts = neko_quad_npts
257 if (this%lgenc) then
258 allocate(htable_i4t2_t::this%facet_map)
259 select type (fmp => this%facet_map)
260 type is(htable_i4t2_t)
261 call fmp%init(this%nelv, facet_data)
262 end select
263
264 allocate(this%facet_neigh(neko_quad_neds, this%nelv))
265
266 call this%hte%init(this%nelv * neko_quad_neds, i)
267 end if
268 else
269 call neko_error("Invalid dimension")
270 end if
271
273 allocate(this%points(this%npts*this%nelv))
274
277 if (this%lgenc) then
278 allocate(this%point_neigh(this%gdim*this%npts*this%nelv))
279 do i = 1, this%gdim*this%npts*this%nelv
280 call this%point_neigh(i)%init()
281 end do
282 end if
283
284 allocate(this%facet_type(2 * this%gdim, this%nelv))
285 this%facet_type = 0
286
287 call this%htp%init(this%npts*this%nelv, i)
288 call this%htel%init(this%nelv, i)
289
290 call this%periodic%init(this%nelv)
291
292 allocate(this%labeled_zones(neko_msh_max_zlbls))
293 do i = 1, neko_msh_max_zlbls
294 call this%labeled_zones(i)%init(this%nelv)
295 end do
296
297 call this%curve%init(this%nelv)
298
299 call distdata_init(this%ddata)
300
301 allocate(this%neigh(0:pe_size-1))
302 this%neigh = .false.
303
304 this%mpts = 0
305 this%mfcs = 0
306 this%meds = 0
307
308 end subroutine mesh_init_common
309
311 subroutine mesh_free(this)
312 class(mesh_t), intent(inout) :: this
313 integer :: i
314
315 call this%htp%free()
316 call this%htf%free()
317 call this%hte%free()
318 call this%htel%free()
319 call distdata_free(this%ddata)
320 call this%curve%free()
321
322 if (allocated(this%dfrmd_el)) then
323 deallocate(this%dfrmd_el)
324 end if
325
326 if (allocated(this%elements)) then
327 do i = 1, this%nelv
328 call this%elements(i)%e%free()
329 deallocate(this%elements(i)%e)
330 end do
331 deallocate(this%elements)
332 end if
333
334 if (allocated(this%facet_map)) then
335 select type (fmp => this%facet_map)
336 type is(htable_i4t2_t)
337 call fmp%free()
338 type is(htable_i4t4_t)
339 call fmp%free()
340 end select
341 deallocate(this%facet_map)
342 end if
343
344 if (allocated(this%facet_neigh)) then
345 deallocate(this%facet_neigh)
346 end if
347
348 if (allocated(this%point_neigh)) then
349 do i = 1, this%gdim * this%npts * this%nelv
350 call this%point_neigh(i)%free()
351 end do
352 deallocate(this%point_neigh)
353 end if
354
355 if (allocated(this%facet_type)) then
356 deallocate(this%facet_type)
357 end if
358 if (allocated(this%labeled_zones)) then
359 do i = 1, neko_msh_max_zlbls
360 call this%labeled_zones(i)%free()
361 end do
362 deallocate(this%labeled_zones)
363 end if
364
365 if (allocated(this%neigh)) then
366 deallocate(this%neigh)
367 end if
368
369 if (allocated(this%neigh_order)) then
370 deallocate(this%neigh_order)
371 end if
372
373 if (allocated(this%points)) then
374 deallocate(this%points)
375 end if
376
377 call this%periodic%free()
378 this%lconn = .false.
379 this%lnumr = .false.
380 this%ldist = .false.
381 this%lgenc = .true.
382
383 end subroutine mesh_free
384
385 subroutine mesh_finalize(this)
386 class(mesh_t), target, intent(inout) :: this
387 integer :: i
388
389 call mesh_generate_flags(this)
390 call mesh_generate_conn(this)
391
392 call this%periodic%finalize()
393 do i = 1, neko_msh_max_zlbls
394 call this%labeled_zones(i)%finalize()
395 end do
396 call this%curve%finalize()
397
398 end subroutine mesh_finalize
399
400 subroutine mesh_generate_flags(this)
401 type(mesh_t), intent(inout) :: this
402 real(kind=dp) :: u(3), v(3), w(3), temp
403 integer :: e
404
405 do e = 1, this%nelv
406 if (this%gdim .eq. 2) then
407 this%dfrmd_el(e) = .false.
408 u = this%elements(e)%e%pts(2)%p%x - this%elements(e)%e%pts(1)%p%x
409 v = this%elements(e)%e%pts(3)%p%x - this%elements(e)%e%pts(1)%p%x
410 temp = u(1)*v(1) + u(2)*v(2)
411 if(.not. abscmp(temp, 0d0)) this%dfrmd_el(e) = .true.
412 else
413 this%dfrmd_el(e) = .false.
414 u = this%elements(e)%e%pts(2)%p%x - this%elements(e)%e%pts(1)%p%x
415 v = this%elements(e)%e%pts(3)%p%x - this%elements(e)%e%pts(1)%p%x
416 w = this%elements(e)%e%pts(5)%p%x - this%elements(e)%e%pts(1)%p%x
417 temp = u(1)*v(1) + u(2)*v(2) + u(3)*v(3)
418 if(.not. abscmp(temp, 0d0)) this%dfrmd_el(e) = .true.
419 temp = u(1)*w(1) + u(2)*w(2) + u(3)*w(3)
420 if(.not. abscmp(temp, 0d0)) this%dfrmd_el(e) = .true.
421 u = this%elements(e)%e%pts(7)%p%x - this%elements(e)%e%pts(8)%p%x
422 v = this%elements(e)%e%pts(6)%p%x - this%elements(e)%e%pts(8)%p%x
423 w = this%elements(e)%e%pts(4)%p%x - this%elements(e)%e%pts(8)%p%x
424 temp = u(1)*v(1) + u(2)*v(2) + u(3)*v(3)
425 if(.not. abscmp(temp, 0d0)) this%dfrmd_el(e) = .true.
426 temp = u(1)*w(1) + u(2)*w(2) + u(3)*w(3)
427 if(.not. abscmp(temp, 0d0)) this%dfrmd_el(e) = .true.
428 end if
429 end do
430 end subroutine mesh_generate_flags
431
433 subroutine mesh_all_deformed(this)
434 class(mesh_t), intent(inout) :: this
435 this%dfrmd_el = .true.
436 end subroutine mesh_all_deformed
437
439 subroutine mesh_generate_conn(this)
440 class(mesh_t), target, intent(inout) :: this
441 type(tuple_i4_t) :: edge
442 type(tuple4_i4_t) :: face, face_comp
443 type(tuple_i4_t) :: facet_data
444 type(stack_i4_t) :: neigh_order
445 class(element_t), pointer :: ep
446 type(tuple_i4_t) :: e
447 type(tuple4_i4_t) :: f
448 integer :: p_local_idx
449 integer :: el, id
450 integer :: i, j, k, ierr, el_glb_idx, n_sides, n_nodes, src, dst
451
452 if (this%lconn) return
453
454 if (.not. this%lgenc) return
455
456 !If we generate connectivity, we do that here.
457 do el = 1, this%nelv
458 ep => this%elements(el)%e
459 select type(ep)
460 type is (hex_t)
461 do i = 1, neko_hex_npts
462 !Only for getting the id
463 call this%add_point(ep%pts(i)%p, id)
464 p_local_idx = this%get_local(this%points(id))
465 !should stack have inout on what we push? would be neat with in
466 id = ep%id()
467 call this%point_neigh(p_local_idx)%push(id)
468 end do
469 do i = 1, neko_hex_nfcs
470 call ep%facet_id(f, i)
471 call this%add_face(f)
472 end do
473
474 do i = 1, neko_hex_neds
475 call ep%edge_id(e, i)
476 call this%add_edge(e)
477 end do
478 type is (quad_t)
479 do i = 1, neko_quad_npts
480 !Only for getting the id
481 call this%add_point(ep%pts(i)%p, id)
482 p_local_idx = this%get_local(this%points(id))
483 !should stack have inout on what we push? would be neat with in
484 id = ep%id()
485 call this%point_neigh(p_local_idx)%push(id)
486 end do
487
488 do i = 1, neko_quad_neds
489 call ep%facet_id(e, i)
490 call this%add_edge(e)
491 end do
492 end select
493 end do
494
495
496 if (this%gdim .eq. 2) then
497 n_sides = 4
498 n_nodes = 2
499 else
500 n_sides = 6
501 n_nodes = 4
502 end if
503
504 ! Compute global number of unique points
505 call mpi_allreduce(this%max_pts_id, this%glb_mpts, 1, &
506 mpi_integer, mpi_max, neko_comm, ierr)
507
508 !
509 ! Find all (local) boundaries
510 !
511
515 select type (fmp => this%facet_map)
516 type is(htable_i4t2_t)
517 do k = 1, 2
518 do i = 1, this%nelv
519 el_glb_idx = i + this%offset_el
520 do j = 1, n_sides
521 call this%elements(i)%e%facet_id(edge, j)
522
523 ! Assume that all facets are on the exterior
524 facet_data%x = [0, 0]
525
526 !check it this face has shown up earlier
527 if (fmp%get(edge, facet_data) .eq. 0) then
528 !if element is already recognized on face
529 if (facet_data%x(1) .eq. el_glb_idx ) then
530 this%facet_neigh(j, i) = facet_data%x(2)
531 else if( facet_data%x(2) .eq. el_glb_idx) then
532 this%facet_neigh(j, i) = facet_data%x(1)
533 !if this is the second element, arrange so low id is first
534 else if(facet_data%x(1) .gt. el_glb_idx) then
535 facet_data%x(2) = facet_data%x(1)
536 facet_data%x(1) = el_glb_idx
537 this%facet_neigh(j, i) = facet_data%x(2)
538 call fmp%set(edge, facet_data)
539 else if(facet_data%x(1) .lt. el_glb_idx) then
540 facet_data%x(2) = el_glb_idx
541 this%facet_neigh(j, i) = facet_data%x(1)
542 call fmp%set(edge, facet_data)
543 end if
544 else
545 facet_data%x(1) = el_glb_idx
546 this%facet_neigh(j, i) = facet_data%x(2)
547 call fmp%set(edge, facet_data)
548 end if
549 end do
550 end do
551 end do
552 type is(htable_i4t4_t)
553
554 do k = 1, 2
555 do i = 1, this%nelv
556 el_glb_idx = i + this%offset_el
557 do j = 1, n_sides
558 call this%elements(i)%e%facet_id(face, j)
559
560 facet_data%x = (/ 0, 0/)
561
562 !check it this face has shown up earlier
563 if (fmp%get(face, facet_data) .eq. 0) then
564 !if element is already recognized on face
565 if (facet_data%x(1) .eq. el_glb_idx ) then
566 this%facet_neigh(j, i) = facet_data%x(2)
567 call this%elements(i)%e%facet_id(face_comp, &
568 j + (2*mod(j, 2) - 1))
569 if (face_comp .eq. face) then
570 facet_data%x(2) = el_glb_idx
571 this%facet_neigh(j, i) = facet_data%x(1)
572 call fmp%set(face, facet_data)
573 end if
574 else if( facet_data%x(2) .eq. el_glb_idx) then
575 this%facet_neigh(j, i) = facet_data%x(1)
576 !if this is the second element, arrange so low id is first
577 else if(facet_data%x(1) .gt. el_glb_idx) then
578 facet_data%x(2) = facet_data%x(1)
579 facet_data%x(1) = el_glb_idx
580 this%facet_neigh(j, i) = facet_data%x(2)
581 call fmp%set(face, facet_data)
582 else if(facet_data%x(1) .lt. el_glb_idx) then
583 facet_data%x(2) = el_glb_idx
584 this%facet_neigh(j, i) = facet_data%x(1)
585 call fmp%set(face, facet_data)
586 end if
587 else
588 facet_data%x(1) = el_glb_idx
589 this%facet_neigh(j, i) = 0
590 call fmp%set(face, facet_data)
591 end if
592 end do
593 end do
594 end do
595 class default
596 call neko_error('Invalid facet map')
597 end select
598
599
600 !
601 ! Find all external (between PEs) boundaries
602 !
603 if (pe_size .gt. 1) then
604
606
607 !
608 ! Generate neighbour exchange order
609 !
610 call neigh_order%init(pe_size)
611
612 do i = 1, pe_size - 1
613 src = modulo(pe_rank - i + pe_size, pe_size)
614 dst = modulo(pe_rank + i, pe_size)
615 if (this%neigh(src) .or. this%neigh(dst)) then
616 j = i ! adhere to standards...
617 call neigh_order%push(j)
618 end if
619 end do
620
621 allocate(this%neigh_order(neigh_order%size()))
622 select type(order => neigh_order%data)
623 type is (integer)
624 do i = 1, neigh_order%size()
625 this%neigh_order(i) = order(i)
626 end do
627 end select
628 call neigh_order%free()
629
631 else
632 allocate(this%neigh_order(1))
633 this%neigh_order = 1
634 end if
635
636 !
637 ! Find all internal/extenral edge connections
638 ! (Note it needs to be called after external point connections has
639 ! been established)
640 !
641 if (this%gdim .eq. 3) then
642 call mesh_generate_edge_conn(this)
643 end if
644
645
647
648 this%lconn = .true.
649
650 end subroutine mesh_generate_conn
651
654 type(mesh_t), intent(inout) :: this
655 type(tuple_i4_t) :: edge, edge2
656 type(tuple4_i4_t) :: face, face2
657 type(tuple_i4_t) :: facet_data
658 type(stack_i4_t) :: buffer
659 type(mpi_status) :: status
660 type(mpi_request) :: send_req, recv_req
661 integer, allocatable :: recv_buffer(:)
662 integer :: i, j, k, el_glb_idx, n_sides, n_nodes, facet, element, l
663 integer :: max_recv, ierr, src, dst, n_recv, recv_side, neigh_el
664
665
666 if (this%gdim .eq. 2) then
667 n_sides = 4
668 n_nodes = 2
669 else
670 n_sides = 6
671 n_nodes = 4
672 end if
673
674 call buffer%init()
675
676 ! Build send buffers containing
677 ! [el_glb_idx, side number, facet_id (global ids of points)]
678 do i = 1, this%nelv
679 el_glb_idx = i + this%offset_el
680 do j = 1, n_sides
681 facet = j ! Adhere to standards...
682 if (this%facet_neigh(j, i) .eq. 0) then
683 if (n_nodes .eq. 2) then
684 call this%elements(i)%e%facet_id(edge, j)
685 call buffer%push(el_glb_idx)
686 call buffer%push(facet)
687 do k = 1, n_nodes
688 call buffer%push(edge%x(k))
689 end do
690 else
691 call this%elements(i)%e%facet_id(face, j)
692 call buffer%push(el_glb_idx)
693 call buffer%push(facet)
694 do k = 1, n_nodes
695 call buffer%push(face%x(k))
696 end do
697 end if
698 end if
699 end do
700 end do
701
702
703 call mpi_allreduce(buffer%size(), max_recv, 1, &
704 mpi_integer, mpi_max, neko_comm, ierr)
705
706 allocate(recv_buffer(max_recv))
707
708 do i = 1, size(this%neigh_order)
709 src = modulo(pe_rank - this%neigh_order(i) + pe_size, pe_size)
710 dst = modulo(pe_rank + this%neigh_order(i), pe_size)
711
712 if (this%neigh(src)) then
713 call mpi_irecv(recv_buffer, max_recv, mpi_integer, &
714 src, 0, neko_comm, recv_req, ierr)
715 end if
716
717 if (this%neigh(dst)) then
718 call mpi_isend(buffer%array(), buffer%size(), mpi_integer, &
719 dst, 0, neko_comm, send_req, ierr)
720 end if
721
722 if (this%neigh(src)) then
723 call mpi_wait(recv_req, status, ierr)
724 call mpi_get_count(status, mpi_integer, n_recv, ierr)
725
726 select type (fmp => this%facet_map)
727 type is(htable_i4t2_t)
728 do j = 1, n_recv, n_nodes + 2
729 neigh_el = recv_buffer(j)
730 recv_side = recv_buffer(j+1)
731
732 edge = (/ recv_buffer(j+2), recv_buffer(j+3) /)
733
734 facet_data = (/ 0, 0 /)
735 !Check if the face is present on this PE
736 if (fmp%get(edge, facet_data) .eq. 0) then
737 element = facet_data%x(1) - this%offset_el
738 !Check which side is connected
739 do l = 1, n_sides
740 call this%elements(element)%e%facet_id(edge2, l)
741 if(edge2 .eq. edge) then
742 facet = l
743 exit
744 end if
745 end do
746 this%facet_neigh(facet, element) = -neigh_el
747 facet_data%x(2) = -neigh_el
748
749 ! Update facet map
750 call fmp%set(edge, facet_data)
751
752 call distdata_set_shared_el_facet(this%ddata, element, facet)
753
754 if (this%hte%get(edge, facet) .eq. 0) then
755 call distdata_set_shared_facet(this%ddata, facet)
756 else
757 call neko_error("Invalid shared edge")
758 end if
759
760 end if
761
762 end do
763 type is(htable_i4t4_t)
764 do j = 1, n_recv, n_nodes + 2
765 neigh_el = recv_buffer(j)
766 recv_side = recv_buffer(j+1)
767
768 face%x = (/ recv_buffer(j+2), recv_buffer(j+3), &
769 recv_buffer(j+4), recv_buffer(j+5) /)
770
771
772 facet_data%x = (/ 0, 0 /)
773
774 !Check if the face is present on this PE
775 if (fmp%get(face, facet_data) .eq. 0) then
776 ! Determine opposite side and update neighbor
777 element = facet_data%x(1) - this%offset_el
778 do l = 1, 6
779 call this%elements(element)%e%facet_id(face2, l)
780 if(face2 .eq. face) then
781 facet = l
782 exit
783 end if
784 end do
785 this%facet_neigh(facet, element) = -neigh_el
786 facet_data%x(2) = -neigh_el
787
788 ! Update facet map
789 call fmp%set(face, facet_data)
790
791 call distdata_set_shared_el_facet(this%ddata, element, facet)
792
793 if (this%htf%get(face, facet) .eq. 0) then
794 call distdata_set_shared_facet(this%ddata, facet)
795 else
796 call neko_error("Invalid shared face")
797 end if
798
799
800 end if
801
802 end do
803 end select
804 end if
805
806 if (this%neigh(dst)) then
807 call mpi_wait(send_req, mpi_status_ignore, ierr)
808 end if
809
810 end do
811
812
813 deallocate(recv_buffer)
814
815 call buffer%free()
816
818
821 type(mesh_t), intent(inout) :: this
822 type(stack_i4_t) :: send_buffer
823 type(mpi_status) :: status
824 integer, allocatable :: recv_buffer(:)
825 integer :: i, j, k
826 integer :: max_recv, ierr, src, dst, n_recv, neigh_el
827 integer :: pt_glb_idx, pt_loc_idx, num_neigh
828 integer, contiguous, pointer :: neighs(:)
829
830
831 call send_buffer%init(this%mpts * 2)
832
833 ! Build send buffers containing
834 ! [pt_glb_idx, #neigh, neigh id_1 ....neigh_id_n]
835 do i = 1, this%mpts
836 pt_glb_idx = this%points(i)%id() ! Adhere to standards...
837 num_neigh = this%point_neigh(i)%size()
838 call send_buffer%push(pt_glb_idx)
839 call send_buffer%push(num_neigh)
840
841 neighs => this%point_neigh(i)%array()
842 do j = 1, this%point_neigh(i)%size()
843 call send_buffer%push(neighs(j))
844 end do
845 end do
846
847 call mpi_allreduce(send_buffer%size(), max_recv, 1, &
848 mpi_integer, mpi_max, neko_comm, ierr)
849 allocate(recv_buffer(max_recv))
850
851 do i = 1, pe_size - 1
852 src = modulo(pe_rank - i + pe_size, pe_size)
853 dst = modulo(pe_rank + i, pe_size)
854
855 call mpi_sendrecv(send_buffer%array(), send_buffer%size(), &
856 mpi_integer, dst, 0, recv_buffer, max_recv, mpi_integer, src, 0, &
857 neko_comm, status, ierr)
858
859 call mpi_get_count(status, mpi_integer, n_recv, ierr)
860
861 j = 1
862 do while (j .le. n_recv)
863 pt_glb_idx = recv_buffer(j)
864 num_neigh = recv_buffer(j + 1)
865 ! Check if the point is present on this PE
866 pt_loc_idx = this%have_point_glb_idx(pt_glb_idx)
867 if (pt_loc_idx .gt. 0) then
868 this%neigh(src) = .true.
869 do k = 1, num_neigh
870 neigh_el = -recv_buffer(j + 1 + k)
871 call this%point_neigh(pt_loc_idx)%push(neigh_el)
872 call distdata_set_shared_point(this%ddata, pt_loc_idx)
873 end do
874 end if
875 j = j + (2 + num_neigh)
876 end do
877
878 end do
879
880 deallocate(recv_buffer)
881 call send_buffer%free()
882
884
888 subroutine mesh_generate_edge_conn(this)
889 type(mesh_t), target, intent(inout) :: this
890 type(htable_iter_i4t2_t) :: it
891 type(tuple_i4_t), pointer :: edge
892 type(uset_i8_t), target :: edge_idx, ghost, owner
893 type(stack_i8_t), target :: send_buff
894 type(htable_i8_t) :: glb_to_loc
895 type(mpi_status) :: status
896 type(mpi_request) :: send_req, recv_req
897 integer, contiguous, pointer :: p1(:), p2(:), ns_id(:)
898 integer :: i, j, id, ierr, num_edge_glb, edge_offset, num_edge_loc
899 integer :: k, l , shared_offset, glb_nshared, n_glb_id
900 integer(kind=i8) :: C, glb_max, glb_id
901 integer(kind=i8), pointer :: glb_ptr
902 integer(kind=i8), allocatable :: recv_buff(:)
903 logical :: shared_edge
904 type(stack_i4_t), target :: non_shared_edges
905 integer :: max_recv, src, dst, n_recv
906
907
909 allocate(this%ddata%local_to_global_edge(this%meds))
910
911 call edge_idx%init(this%hte%num_entries())
912 call send_buff%init(this%hte%num_entries())
913 call owner%init(this%hte%num_entries())
914
915 call glb_to_loc%init(32, i)
916
917 !
918 ! Determine/ constants used to generate unique global edge numbers
919 ! for shared edges
920 !
921 c = int(this%glb_nelv, i8) * int(neko_hex_neds, i8)
922
923 num_edge_glb = 2* this%meds
924 call mpi_allreduce(mpi_in_place, num_edge_glb, 1, &
925 mpi_integer, mpi_sum, neko_comm, ierr)
926
927 glb_max = int(num_edge_glb, i8)
928
929 call non_shared_edges%init(this%hte%num_entries())
930
931 call it%init(this%hte)
932 do while(it%next())
933 edge => it%key()
934 call it%data(id)
935
936 k = this%have_point_glb_idx(edge%x(1))
937 l = this%have_point_glb_idx(edge%x(2))
938 p1 => this%point_neigh(k)%array()
939 p2 => this%point_neigh(l)%array()
940
941 shared_edge = .false.
942
943 ! Find edge neighbor from point neighbors
944 do i = 1, this%point_neigh(k)%size()
945 do j = 1, this%point_neigh(l)%size()
946 if ((p1(i) .eq. p2(j)) .and. &
947 (p1(i) .lt. 0) .and. (p2(j) .lt. 0)) then
948 call distdata_set_shared_edge(this%ddata, id)
949 shared_edge = .true.
950 end if
951 end do
952 end do
953
954 ! Generate a unique id for the shared edge as,
955 ! ((e1 * C) + e2 )) + glb_max if e1 > e2
956 ! ((e2 * C) + e1 )) + glb_max if e2 > e1
957 if (shared_edge) then
958 glb_id = ((int(edge%x(1), i8)) + int(edge%x(2), i8)*c) + glb_max
959 call glb_to_loc%set(glb_id, id)
960 call edge_idx%add(glb_id)
961 call owner%add(glb_id) ! Always assume the PE is the owner
962 call send_buff%push(glb_id)
963 else
964 call non_shared_edges%push(id)
965 end if
966 end do
967
968 ! Determine start offset for global numbering of locally owned edges
969 edge_offset = 0
970 num_edge_loc = non_shared_edges%size()
971 call mpi_exscan(num_edge_loc, edge_offset, 1, &
972 mpi_integer, mpi_sum, neko_comm, ierr)
973 edge_offset = edge_offset + 1
974
975 ! Construct global numbering of locally owned edges
976 ns_id => non_shared_edges%array()
977 do i = 1, non_shared_edges%size()
978 call distdata_set_local_to_global_edge(this%ddata, ns_id(i), edge_offset)
979 edge_offset = edge_offset + 1
980 end do
981 nullify(ns_id)
982
983 !
984 ! Renumber shared edges into integer range
985 !
986
987 call mpi_allreduce(send_buff%size(), max_recv, 1, &
988 mpi_integer, mpi_max, neko_comm, ierr)
989
990 call ghost%init(send_buff%size())
991
992 allocate(recv_buff(max_recv))
993
994 do i = 1, size(this%neigh_order)
995 src = modulo(pe_rank - this%neigh_order(i) + pe_size, pe_size)
996 dst = modulo(pe_rank + this%neigh_order(i), pe_size)
997
998 if (this%neigh(src)) then
999 call mpi_irecv(recv_buff, max_recv, mpi_integer8, &
1000 src, 0, neko_comm, recv_req, ierr)
1001 end if
1002
1003 if (this%neigh(dst)) then
1004 ! We should use the %array() procedure, which works great for
1005 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
1006 ! certain data types
1007 select type(sbarray=>send_buff%data)
1008 type is (integer(i8))
1009 call mpi_isend(sbarray, send_buff%size(), mpi_integer8, &
1010 dst, 0, neko_comm, send_req, ierr)
1011 end select
1012 end if
1013
1014 if (this%neigh(src)) then
1015 call mpi_wait(recv_req, status, ierr)
1016 call mpi_get_count(status, mpi_integer8, n_recv, ierr)
1017
1018 do j = 1, n_recv
1019 if ((edge_idx%element(recv_buff(j))) .and. (src .lt. pe_rank)) then
1020 call ghost%add(recv_buff(j))
1021 call owner%remove(recv_buff(j))
1022 end if
1023 end do
1024 end if
1025
1026 if (this%neigh(dst)) then
1027 call mpi_wait(send_req, mpi_status_ignore, ierr)
1028 end if
1029 end do
1030
1031
1032 ! Determine start offset for global numbering of shared edges
1033 glb_nshared = num_edge_loc
1034 call mpi_allreduce(mpi_in_place, glb_nshared, 1, &
1035 mpi_integer, mpi_sum, neko_comm, ierr)
1036
1037 shared_offset = 0
1038 call mpi_exscan(owner%size(), shared_offset, 1, &
1039 mpi_integer, mpi_sum, neko_comm, ierr)
1040 shared_offset = shared_offset + glb_nshared + 1
1041
1042 ! Renumber locally owned set of shared edges
1043 call send_buff%clear()
1044 call owner%iter_init()
1045 do while (owner%iter_next())
1046 glb_ptr => owner%iter_value()
1047 if (glb_to_loc%get(glb_ptr, id) .eq. 0) then
1048 call distdata_set_local_to_global_edge(this%ddata, id, shared_offset)
1049
1050 ! Add new number to send data as [old_glb_id new_glb_id] for each edge
1051 call send_buff%push(glb_ptr) ! Old glb_id integer*8
1052 glb_id = int(shared_offset, i8) ! Waste some space here...
1053 call send_buff%push(glb_id) ! New glb_id integer*4
1054
1055 shared_offset = shared_offset + 1
1056 else
1057 call neko_error('Invalid edge id')
1058 end if
1059 end do
1060 nullify(glb_ptr)
1061
1062 ! Determine total number of unique edges in the mesh
1063 ! (This can probably be done in a clever way...)
1064 this%glb_meds = shared_offset -1
1065 call mpi_allreduce(mpi_in_place, this%glb_meds, 1, &
1066 mpi_integer, mpi_max, neko_comm, ierr)
1067
1068 !
1069 ! Update ghosted edges with new global id
1070 !
1071
1072 call mpi_allreduce(send_buff%size(), max_recv, 1, &
1073 mpi_integer, mpi_max, neko_comm, ierr)
1074
1075 deallocate(recv_buff)
1076 allocate(recv_buff(max_recv))
1077
1078
1079 do i = 1, size(this%neigh_order)
1080 src = modulo(pe_rank - this%neigh_order(i) + pe_size, pe_size)
1081 dst = modulo(pe_rank + this%neigh_order(i), pe_size)
1082
1083 if (this%neigh(src)) then
1084 call mpi_irecv(recv_buff, max_recv, mpi_integer8, &
1085 src, 0, neko_comm, recv_req, ierr)
1086 end if
1087
1088 if (this%neigh(dst)) then
1089 ! We should use the %array() procedure, which works great for
1090 ! GNU, Intel and NEC, but it breaks horribly on Cray when using
1091 ! certain data types
1092 select type(sbarray=>send_buff%data)
1093 type is (integer(i8))
1094 call mpi_isend(sbarray, send_buff%size(), mpi_integer8, &
1095 dst, 0, neko_comm, send_req, ierr)
1096 end select
1097 end if
1098
1099 if (this%neigh(src)) then
1100 call mpi_wait(recv_req, status, ierr)
1101 call mpi_get_count(status, mpi_integer8, n_recv, ierr)
1102
1103 do j = 1, n_recv, 2
1104 if (ghost%element(recv_buff(j))) then
1105 if (glb_to_loc%get(recv_buff(j), id) .eq. 0) then
1106 n_glb_id = int(recv_buff(j + 1 ), 4)
1107 call distdata_set_local_to_global_edge(this%ddata, id, &
1108 n_glb_id)
1109 else
1110 call neko_error('Invalid edge id')
1111 end if
1112 end if
1113 end do
1114 end if
1115
1116 if (this%neigh(dst)) then
1117 call mpi_wait(send_req, mpi_status_ignore, ierr)
1118 end if
1119 end do
1120
1121 deallocate(recv_buff)
1122 call glb_to_loc%free()
1123 call send_buff%free()
1124 call edge_idx%free()
1125 call non_shared_edges%free()
1126 call ghost%free()
1127 call owner%free()
1128
1129 end subroutine mesh_generate_edge_conn
1130
1133 type(mesh_t), target, intent(inout) :: this
1134 type(htable_iter_i4t4_t), target :: face_it
1135 type(htable_iter_i4t2_t), target :: edge_it
1136 type(tuple4_i4_t), pointer :: face, fd(:)
1137 type(tuple_i4_t), pointer :: edge, ed(:)
1138 type(tuple_i4_t) :: facet_data
1139 type(tuple4_i4_t) :: recv_face
1140 type(tuple_i4_t) :: recv_edge
1141 type(stack_i4t4_t) :: face_owner
1142 type(htable_i4t4_t) :: face_ghost
1143 type(stack_i4t2_t) :: edge_owner
1144 type(htable_i4t2_t) :: edge_ghost
1145 type(stack_i4_t) :: send_buff
1146 type(mpi_status) :: status
1147 type(mpi_request) :: send_req, recv_req
1148 integer, allocatable :: recv_buff(:)
1149 integer :: non_shared_facets, shared_facets, facet_offset
1150 integer :: id, glb_nshared, shared_offset, owned_facets
1151 integer :: i, j, ierr, max_recv, src, dst, n_recv
1152
1153 shared_facets = this%ddata%shared_facet%size()
1154
1156 if (this%gdim .eq. 2) then
1157 allocate(this%ddata%local_to_global_facet(this%meds))
1158 call edge_owner%init(this%meds)
1159 call edge_ghost%init(64, i)
1160 non_shared_facets = this%hte%num_entries() - shared_facets
1161 else
1162 allocate(this%ddata%local_to_global_facet(this%mfcs))
1163 call face_owner%init(this%mfcs)
1164 call face_ghost%init(64, i)
1165 non_shared_facets = this%htf%num_entries() - shared_facets
1166 end if
1167
1169
1170 facet_offset = 0
1171 call mpi_exscan(non_shared_facets, facet_offset, 1, &
1172 mpi_integer, mpi_sum, neko_comm, ierr)
1173 facet_offset = facet_offset + 1
1174
1175 ! Determine ownership of shared facets
1176 if (this%gdim .eq. 2) then
1177 call edge_it%init(this%hte)
1178 do while (edge_it%next())
1179 call edge_it%data(id)
1180 edge => edge_it%key()
1181 if (.not. this%ddata%shared_facet%element(id)) then
1182 call distdata_set_local_to_global_facet(this%ddata, &
1183 id, facet_offset)
1184 facet_offset = facet_offset + 1
1185 else
1186 select type(fmp => this%facet_map)
1187 type is(htable_i4t2_t)
1188 if (fmp%get(edge, facet_data) .eq. 0) then
1189 if (facet_data%x(2) .lt. 0) then
1190 if (abs(facet_data%x(2)) .lt. (this%offset_el + 1)) then
1191 call edge_ghost%set(edge, id)
1192 else
1193 call edge_owner%push(edge)
1194 end if
1195 else
1196 call neko_error("Invalid edge neigh.")
1197 end if
1198 end if
1199 end select
1200 end if
1201 end do
1202 owned_facets = edge_owner%size()
1203 else
1204 call face_it%init(this%htf)
1205 do while (face_it%next())
1206 call face_it%data(id)
1207 face => face_it%key()
1208 if (.not. this%ddata%shared_facet%element(id)) then
1209 call distdata_set_local_to_global_facet(this%ddata, &
1210 id, facet_offset)
1211 facet_offset = facet_offset + 1
1212 else
1213 select type(fmp => this%facet_map)
1214 type is(htable_i4t4_t)
1215 if (fmp%get(face, facet_data) .eq. 0) then
1216 if (facet_data%x(2) .lt. 0) then
1217 if (abs(facet_data%x(2)) .lt. (this%offset_el + 1)) then
1218 call face_ghost%set(face, id)
1219 else
1220 call face_owner%push(face)
1221 end if
1222 else
1223 call neko_error("Invalid face neigh.")
1224 end if
1225 end if
1226 end select
1227 end if
1228 end do
1229 owned_facets = face_owner%size()
1230 end if
1231
1232 ! Determine start offset for global numbering of shared facets
1233 glb_nshared = non_shared_facets
1234 call mpi_allreduce(mpi_in_place, glb_nshared, 1, &
1235 mpi_integer, mpi_sum, neko_comm, ierr)
1236
1237 shared_offset = 0
1238 call mpi_exscan(owned_facets, shared_offset, 1, &
1239 mpi_integer, mpi_sum, neko_comm, ierr)
1240 shared_offset = shared_offset + glb_nshared + 1
1241
1242 if (this%gdim .eq. 2) then
1243
1244 if (owned_facets .gt. 32) then
1245 call send_buff%init(owned_facets)
1246 else
1247 call send_buff%init()
1248 end if
1249
1250 ed => edge_owner%array()
1251 do i = 1, edge_owner%size()
1252 if (this%hte%get(ed(i), id) .eq. 0) then
1253 call distdata_set_local_to_global_facet(this%ddata, id, &
1254 shared_offset)
1255
1256 ! Add new number to send buffer
1257 ! [edge id1 ... edge idn new_glb_id]
1258 do j = 1, 2
1259 call send_buff%push(ed(i)%x(j))
1260 end do
1261 call send_buff%push(shared_offset)
1262
1263 shared_offset = shared_offset + 1
1264 end if
1265 end do
1266
1267 else
1268
1269 if (owned_facets .gt. 32) then
1270 call send_buff%init(owned_facets)
1271 else
1272 call send_buff%init()
1273 end if
1274
1275 fd => face_owner%array()
1276 do i = 1, face_owner%size()
1277 if (this%htf%get(fd(i), id) .eq. 0) then
1278 call distdata_set_local_to_global_facet(this%ddata, id, &
1279 shared_offset)
1280
1281 ! Add new number to send buffer
1282 ! [face id1 ... face idn new_glb_id]
1283 do j = 1, 4
1284 call send_buff%push(fd(i)%x(j))
1285 end do
1286 call send_buff%push(shared_offset)
1287
1288 shared_offset = shared_offset + 1
1289 end if
1290 end do
1291 nullify(fd)
1292
1293 end if
1294
1295 ! Determine total number of unique facets in the mesh
1296 ! (This can probably be done in a clever way...)
1297 this%glb_mfcs = shared_offset - 1
1298 call mpi_allreduce(mpi_in_place, this%glb_mfcs, 1, &
1299 mpi_integer, mpi_max, neko_comm, ierr)
1300
1301 !
1302 ! Update ghosted facets with new global id
1303 !
1304
1305 call mpi_allreduce(send_buff%size(), max_recv, 1, &
1306 mpi_integer, mpi_max, neko_comm, ierr)
1307
1308 allocate(recv_buff(max_recv))
1309
1311 do i = 1, size(this%neigh_order)
1312 src = modulo(pe_rank - this%neigh_order(i) + pe_size, pe_size)
1313 dst = modulo(pe_rank + this%neigh_order(i), pe_size)
1314
1315 if (this%neigh(src)) then
1316 call mpi_irecv(recv_buff, max_recv, mpi_integer, &
1317 src, 0, neko_comm, recv_req, ierr)
1318 end if
1319
1320 if (this%neigh(dst)) then
1321 call mpi_isend(send_buff%array(), send_buff%size(), mpi_integer, &
1322 dst, 0, neko_comm, send_req, ierr)
1323 end if
1324
1325 if (this%neigh(src)) then
1326 call mpi_wait(recv_req, status, ierr)
1327 call mpi_get_count(status, mpi_integer, n_recv, ierr)
1328
1329 if (this%gdim .eq. 2) then
1330 do j = 1, n_recv, 3
1331
1332 recv_edge = (/recv_buff(j), recv_buff(j+1)/)
1333
1334 ! Check if the PE has the shared edge
1335 if (edge_ghost%get(recv_edge, id) .eq. 0) then
1336 call distdata_set_local_to_global_facet(this%ddata, &
1337 id, recv_buff(j+2))
1338 end if
1339 end do
1340 else
1341 do j = 1, n_recv, 5
1342
1343 recv_face = (/recv_buff(j), recv_buff(j+1), &
1344 recv_buff(j+2), recv_buff(j+3) /)
1345
1346 ! Check if the PE has the shared face
1347 if (face_ghost%get(recv_face, id) .eq. 0) then
1348 call distdata_set_local_to_global_facet(this%ddata, &
1349 id, recv_buff(j+4))
1350 end if
1351 end do
1352 end if
1353 end if
1354
1355 if (this%neigh(dst)) then
1356 call mpi_wait(send_req, mpi_status_ignore, ierr)
1357 end if
1358
1359 end do
1360
1361 if (this%gdim .eq. 2) then
1362 call edge_owner%free()
1363 call edge_ghost%free()
1364 else
1365 call face_owner%free()
1366 call face_ghost%free()
1367 end if
1368
1369 call send_buff%free()
1370 deallocate(recv_buff)
1371
1372 end subroutine mesh_generate_facet_numbering
1373
1374
1376 subroutine mesh_add_quad(this, el, el_glb, p1, p2, p3, p4)
1377 class(mesh_t), target, intent(inout) :: this
1378 integer, value :: el, el_glb
1379 type(point_t), target, intent(inout) :: p1, p2, p3, p4
1380 integer :: p(4)
1381 type(tuple_i4_t) :: e
1382
1383 ! Connectivity invalidated if a new element is added
1384 this%lconn = .false.
1385
1386 ! Numbering invalidated if a new element is added
1387 this%lnumr = .false.
1388
1389 call this%add_point(p1, p(1))
1390 call this%add_point(p2, p(2))
1391 call this%add_point(p3, p(3))
1392 call this%add_point(p4, p(4))
1393
1394 select type (ep => this%elements(el)%e)
1395 type is (quad_t)
1396 call ep%init(el_glb, &
1397 this%points(p(1)), this%points(p(2)), &
1398 this%points(p(3)), this%points(p(4)))
1399
1400
1401 class default
1402 call neko_error('Invalid element type')
1403 end select
1404
1405 end subroutine mesh_add_quad
1406
1408 subroutine mesh_add_hex(this, el, el_glb, p1, p2, p3, p4, p5, p6, p7, p8)
1409 class(mesh_t), target, intent(inout) :: this
1410 integer, value :: el, el_glb
1411 type(point_t), target, intent(inout) :: p1, p2, p3, p4, p5, p6, p7, p8
1412 integer :: p(8)
1413 type(tuple4_i4_t) :: f
1414 type(tuple_i4_t) :: e
1415
1416 ! Connectivity invalidated if a new element is added
1417 this%lconn = .false.
1418
1419 ! Numbering invalidated if a new element is added
1420 this%lnumr = .false.
1421
1422 call this%add_point(p1, p(1))
1423 call this%add_point(p2, p(2))
1424 call this%add_point(p3, p(3))
1425 call this%add_point(p4, p(4))
1426 call this%add_point(p5, p(5))
1427 call this%add_point(p6, p(6))
1428 call this%add_point(p7, p(7))
1429 call this%add_point(p8, p(8))
1430
1431 ! Global to local mapping
1432 call this%htel%set(el_glb, el)
1433
1434 select type (ep => this%elements(el)%e)
1435 type is (hex_t)
1436 call ep%init(el_glb, &
1437 this%points(p(1)), this%points(p(2)), &
1438 this%points(p(3)), this%points(p(4)), &
1439 this%points(p(5)), this%points(p(6)), &
1440 this%points(p(7)), this%points(p(8)))
1441 class default
1442 call neko_error('Invalid element type')
1443 end select
1444
1445 end subroutine mesh_add_hex
1446
1448 subroutine mesh_add_point(this, p, idx)
1449 class(mesh_t), intent(inout) :: this
1450 type(point_t), intent(inout) :: p
1451 integer, intent(inout) :: idx
1452 integer :: tmp
1453
1454 tmp = p%id()
1455
1456 this%max_pts_id = max(this%max_pts_id, tmp)
1457
1458 if (tmp .le. 0) then
1459 call neko_error("Invalid point id")
1460 end if
1461
1462 if (this%htp%get(tmp, idx) .gt. 0) then
1463 this%mpts = this%mpts + 1
1464 call this%htp%set(tmp, this%mpts)
1465 this%points(this%mpts) = p
1466 idx = this%mpts
1467 end if
1468
1469 end subroutine mesh_add_point
1470
1472 subroutine mesh_add_face(this, f)
1473 class(mesh_t), intent(inout) :: this
1474 type(tuple4_i4_t), intent(inout) :: f
1475 integer :: idx
1476
1477 if (this%htf%get(f, idx) .gt. 0) then
1478 this%mfcs = this%mfcs + 1
1479 call this%htf%set(f, this%mfcs)
1480 end if
1481
1482 end subroutine mesh_add_face
1483
1485 subroutine mesh_add_edge(this, e)
1486 class(mesh_t), intent(inout) :: this
1487 type(tuple_i4_t), intent(inout) :: e
1488 integer :: idx
1489
1490 if (this%hte%get(e, idx) .gt. 0) then
1491 this%meds = this%meds + 1
1492 call this%hte%set(e, this%meds)
1493 end if
1494
1495 end subroutine mesh_add_edge
1496
1498 subroutine mesh_mark_curve_element(this, e, curve_data, curve_type)
1499 class(mesh_t), intent(inout) :: this
1500 integer, intent(in) :: e
1501 real(kind=dp), dimension(5,12), intent(in) :: curve_data
1502 integer, dimension(12), intent(in) :: curve_type
1503
1504 if (e .gt. this%nelv) then
1505 call neko_error('Invalid element index')
1506 end if
1507 if ((this%gdim .eq. 2 .and. sum(curve_type(5:8)) .gt. 0) ) then
1508 call neko_error('Invalid curve element')
1509 end if
1510 call this%curve%add_element(e, curve_data, curve_type)
1511
1512 end subroutine mesh_mark_curve_element
1513
1515 subroutine mesh_mark_labeled_facet(this, f, e, label)
1516 class(mesh_t), intent(inout) :: this
1517 integer, intent(in) :: f
1518 integer, intent(in) :: e
1519 integer, intent(in) :: label
1520
1521 if (e .gt. this%nelv) then
1522 call neko_error('Invalid element index')
1523 end if
1524
1525 if ((this%gdim .eq. 2 .and. f .gt. 4) .or. &
1526 (this%gdim .eq. 3 .and. f .gt. 6)) then
1527 call neko_error('Invalid facet index')
1528 end if
1529 call this%labeled_zones(label)%add_facet(f, e)
1530 this%facet_type(f,e) = -label
1531
1532 end subroutine mesh_mark_labeled_facet
1533
1535 subroutine mesh_mark_periodic_facet(this, f, e, pf, pe, pids)
1536 class(mesh_t), intent(inout) :: this
1537 integer, intent(in) :: f
1538 integer, intent(in) :: e
1539 integer, intent(in) :: pf
1540 integer, intent(in) :: pe
1541 integer, intent(inout) :: pids(4)
1542 integer, dimension(4) :: org_ids
1543
1544 call this%get_facet_ids(f, e, org_ids)
1545 call this%periodic%add_periodic_facet(f, e, pf, pe, pids, org_ids)
1546 end subroutine mesh_mark_periodic_facet
1547
1549 subroutine mesh_get_facet_ids(this, f, e, pids)
1550 class(mesh_t), intent(inout) :: this
1551 integer, intent(in) :: f
1552 integer, intent(in) :: e
1553 integer, intent(inout) :: pids(4)
1554 type(point_t), pointer :: pi
1555 type(tuple4_i4_t) :: t
1556 type(tuple_i4_t) :: t2
1557
1558 select type(ele => this%elements(e)%e)
1559 type is(hex_t)
1560 call ele%facet_order(t,f)
1561 pids = t%x
1562 type is(quad_t)
1563 call ele%facet_order(t2,f)
1564 pids(1) = t2%x(1)
1565 pids(2) = t2%x(2)
1566 pids(3) = 0
1567 pids(4) = 0
1568 end select
1569 end subroutine mesh_get_facet_ids
1570
1573 class(mesh_t), intent(inout) :: this
1574 integer :: i,j
1575 integer :: f
1576 integer :: e
1577 integer :: pf
1578 integer :: pe
1579 integer :: org_ids(4), pids(4)
1580 type(point_t), pointer :: pi
1581 integer, dimension(4, 6) :: face_nodes = reshape([ &
1582 1,5,7,3, &
1583 2,6,8,4, &
1584 1,2,6,5, &
1585 3,4,8,7, &
1586 1,2,4,3, &
1587 5,6,8,7],&
1588 [4,6])
1589 integer, dimension(2, 4) :: edge_nodes = reshape([ &
1590 1,3, &
1591 2,4, &
1592 1,2, &
1593 3,4],&
1594 [2,4])
1595
1596 do i = 1, this%periodic%size
1597 e = this%periodic%facet_el(i)%x(2)
1598 f = this%periodic%facet_el(i)%x(1)
1599 pe = this%periodic%p_facet_el(i)%x(2)
1600 pf = this%periodic%p_facet_el(i)%x(1)
1601 pids = this%periodic%p_ids(i)%x
1602 call this%get_facet_ids(f, e, pids)
1603 this%periodic%p_ids(i)%x = pids
1604 end do
1605 do i = 1, this%periodic%size
1606 e = this%periodic%facet_el(i)%x(2)
1607 f = this%periodic%facet_el(i)%x(1)
1608 org_ids = this%periodic%org_ids(i)%x
1609 select type(ele => this%elements(e)%e)
1610 type is(hex_t)
1611 do j = 1, 4
1612 pi => ele%pts(face_nodes(j,f))%p
1613 call pi%set_id(org_ids(j))
1614 end do
1615 type is(quad_t)
1616 do j = 1, 2
1617 pi => ele%pts(edge_nodes(j,f))%p
1618 call pi%set_id(org_ids(j))
1619 end do
1620 end select
1621 end do
1622 end subroutine mesh_reset_periodic_ids
1623
1625 subroutine mesh_create_periodic_ids(this, f, e, pf, pe)
1626 class(mesh_t), intent(inout) :: this
1627 integer, intent(in) :: f
1628 integer, intent(in) :: e
1629 integer, intent(in) :: pf
1630 integer, intent(in) :: pe
1631 type(point_t), pointer :: pi, pj
1632 real(kind=dp) :: l(3)
1633 integer :: i, j, id, p_local_idx, match
1634 type(tuple4_i4_t) :: ft
1635 type(tuple_i4_t) :: et
1636 integer, dimension(4, 6) :: face_nodes = reshape([&
1637 1,5,7,3,&
1638 2,6,8,4,&
1639 1,2,6,5,&
1640 3,4,8,7,&
1641 1,2,4,3,&
1642 5,6,8,7],&
1643 [4,6])
1644 integer, dimension(2, 4) :: edge_nodes = reshape([&
1645 1,3,&
1646 2,4,&
1647 1,2,&
1648 3,4 ],&
1649 [2,4])
1650
1651 select type(ele => this%elements(e)%e)
1652 type is(hex_t)
1653 select type(elp => this%elements(pe)%e)
1654 type is(hex_t)
1655 l = 0d0
1656 do i = 1, 4
1657 l = l + ele%pts(face_nodes(i,f))%p%x(1:3) - &
1658 elp%pts(face_nodes(i,pf))%p%x(1:3)
1659 end do
1660 l = l/4
1661 do i = 1, 4
1662 pi => ele%pts(face_nodes(i,f))%p
1663 match = 0
1664 do j = 1, 4
1665 pj => elp%pts(face_nodes(j,pf))%p
1666 if (norm2(pi%x(1:3) - pj%x(1:3) - l) .lt. 1d-7) then
1667 id = min(pi%id(), pj%id())
1668 call pi%set_id(id)
1669 call pj%set_id(id)
1670 p_local_idx = this%get_local(this%points(id))
1671 match = match + 1
1672 end if
1673 end do
1674 if ( match .gt. 1) then
1675 call neko_error('Multiple matches when creating periodic ids')
1676 else if (match .eq. 0) then
1677 call neko_error('Cannot find matching periodic point')
1678 end if
1679 end do
1680 end select
1681 type is(quad_t)
1682 select type(elp => this%elements(pe)%e)
1683 type is(quad_t)
1684 l = 0d0
1685 do i = 1, 2
1686 l = l + ele%pts(edge_nodes(i,f))%p%x(1:3) - &
1687 elp%pts(edge_nodes(i,pf))%p%x(1:3)
1688 end do
1689 l = l/2
1690 do i = 1, 2
1691 pi => ele%pts(edge_nodes(i,f))%p
1692 do j = 1, 2
1693 pj => elp%pts(edge_nodes(j,pf))%p
1694 !whatabout thie tolerance?
1695 if (norm2(pi%x(1:3) - pj%x(1:3) - l) .lt. 1d-7) then
1696 id = min(pi%id(), pj%id())
1697 call pi%set_id(id)
1698 call pj%set_id(id)
1699 p_local_idx = this%get_local(this%points(id))
1700 end if
1701 end do
1702 end do
1703 end select
1704 end select
1705 end subroutine mesh_create_periodic_ids
1706
1709 subroutine mesh_apply_periodic_facet(this, f, e, pf, pe, pids)
1710 class(mesh_t), intent(inout) :: this
1711 integer, intent(in) :: f
1712 integer, intent(in) :: e
1713 integer, intent(in) :: pf
1714 integer, intent(in) :: pe
1715 integer, intent(inout) :: pids(4)
1716 type(point_t), pointer :: pi
1717 integer :: i, id, p_local_idx
1718 type(tuple4_i4_t) :: ft
1719 type(tuple_i4_t) :: et
1720 integer, dimension(4, 6) :: face_nodes = reshape([&
1721 1,5,7,3,&
1722 2,6,8,4,&
1723 1,2,6,5,&
1724 3,4,8,7,&
1725 1,2,4,3,&
1726 5,6,8,7],&
1727 [4,6])
1728 select type(ele => this%elements(e)%e)
1729 type is(hex_t)
1730 do i = 1, 4
1731 pi => ele%pts(face_nodes(i,f))%p
1732 call pi%set_id(pids(i))
1733 call this%add_point(pi, id)
1734 p_local_idx = this%get_local(this%points(id))
1735 end do
1736 end select
1737
1738 end subroutine mesh_apply_periodic_facet
1739
1741 function mesh_get_local_point(this, p) result(local_id)
1742 class(mesh_t), intent(inout) :: this
1743 type(point_t), intent(inout) :: p
1744 integer :: local_id
1745 integer :: tmp
1746
1748 tmp = p%id()
1749
1750 if (this%htp%get(tmp, local_id) .gt. 0) then
1751 call neko_error('Invalid global id (local point)')
1752 end if
1753
1754 end function mesh_get_local_point
1755
1758 function mesh_get_local_edge(this, e) result(local_id)
1759 class(mesh_t), intent(inout) :: this
1760 type(tuple_i4_t), intent(inout) :: e
1761 integer :: local_id
1762
1763 if (this%hte%get(e, local_id) .gt. 0) then
1764 call neko_error('Invalid global id (local edge)')
1765 end if
1766
1767 end function mesh_get_local_edge
1768
1770 function mesh_get_local_facet(this, f) result(local_id)
1771 class(mesh_t), intent(inout) :: this
1772 type(tuple4_i4_t), intent(inout) :: f
1773 integer :: local_id
1774
1775 if (this%htf%get(f, local_id) .gt. 0) then
1776 call neko_error('Invalid global id (local facet)')
1777 end if
1778
1779 end function mesh_get_local_facet
1780
1782 function mesh_get_global_edge(this, e) result(global_id)
1783 class(mesh_t), intent(inout) :: this
1784 type(tuple_i4_t), intent(inout) :: e
1785 integer :: global_id
1786
1787 global_id = this%get_local(e)
1788
1789 if (this%gdim .eq. 2) then
1790 if (pe_size .gt. 1) then
1791 global_id = this%ddata%local_to_global_facet(global_id)
1792 end if
1793 else
1794 if (pe_size .gt. 1) then
1795 global_id = this%ddata%local_to_global_edge(global_id)
1796 end if
1797 end if
1798
1799 end function mesh_get_global_edge
1800
1802 function mesh_get_global_facet(this, f) result(global_id)
1803 class(mesh_t), intent(inout) :: this
1804 type(tuple4_i4_t), intent(inout) :: f
1805 integer :: global_id
1806
1807 global_id = this%get_local_facet(f)
1808
1809 if (pe_size .gt. 1) then
1810 global_id = this%ddata%local_to_global_facet(global_id)
1811 end if
1812
1813 end function mesh_get_global_facet
1814
1815
1819 function mesh_have_point_glb_idx(this, index) result(local_id)
1820 class(mesh_t), intent(inout) :: this
1821 integer, intent(inout) :: index
1822 integer :: local_id
1823
1824 if (this%htp%get(index, local_id) .eq. 1) then
1825 local_id = -1
1826 end if
1827
1828 end function mesh_have_point_glb_idx
1829
1830
1832 function mesh_is_shared_point(this, p) result(shared)
1833 class(mesh_t), intent(inout) :: this
1834 type(point_t), intent(inout) :: p
1835 integer :: local_index
1836 logical shared
1837
1838 local_index = this%get_local(p)
1839 shared = this%ddata%shared_point%element(local_index)
1840
1841 end function mesh_is_shared_point
1842
1843
1846 function mesh_is_shared_edge(this, e) result(shared)
1847 class(mesh_t), intent(inout) :: this
1848 type(tuple_i4_t), intent(inout) :: e
1849 integer :: local_index
1850 logical shared
1851 local_index = this%get_local(e)
1852 if (this%gdim .eq. 2) then
1853 shared = this%ddata%shared_facet%element(local_index)
1854 else
1855 shared = this%ddata%shared_edge%element(local_index)
1856 end if
1857 end function mesh_is_shared_edge
1858
1860 function mesh_is_shared_facet(this, f) result(shared)
1861 class(mesh_t), intent(inout) :: this
1862 type(tuple4_i4_t), intent(inout) :: f
1863 integer :: local_index
1864 logical shared
1865
1866 local_index = this%get_local(f)
1867 shared = this%ddata%shared_facet%element(local_index)
1868
1869 end function mesh_is_shared_facet
1870
1871end module mesh
Definition comm.F90:1
integer pe_rank
MPI rank.
Definition comm.F90:51
type(mpi_comm) neko_comm
MPI communicator.
Definition comm.F90:38
integer pe_size
MPI size of communicator.
Definition comm.F90:54
Defines a domain as a subset of facets in a mesh.
Definition curve.f90:2
Defines practical data distributions.
Definition datadist.f90:34
Distributed mesh data.
Definition distdata.f90:34
subroutine, public distdata_set_shared_point(ddata, point)
Mark a point as shared.
Definition distdata.f90:127
subroutine, public distdata_set_shared_el_facet(ddata, element, side)
Mark an element's facet as shared.
Definition distdata.f90:96
subroutine, public distdata_set_local_to_global_facet(ddata, local, global)
Set local to global mapping (facets)
Definition distdata.f90:136
subroutine, public distdata_free(ddata)
Free a distdata type.
Definition distdata.f90:76
subroutine, public distdata_set_local_to_global_edge(ddata, local, global)
Set local to global mapping (edges)
Definition distdata.f90:146
subroutine, public distdata_init(ddata)
Initialise a distdata type.
Definition distdata.f90:62
subroutine, public distdata_set_shared_facet(ddata, facet)
Mark a facet as shared.
Definition distdata.f90:108
subroutine, public distdata_set_shared_edge(ddata, edge)
Mark an element's edge as shared.
Definition distdata.f90:118
Defines a zone as a subset of facets in a mesh.
Defines a hexahedron element.
Definition hex.f90:34
integer, parameter, public neko_hex_npts
Number of points.
Definition hex.f90:42
integer, parameter, public neko_hex_nfcs
Number of faces.
Definition hex.f90:43
integer, parameter, public neko_hex_neds
Number of edges.
Definition hex.f90:44
Implements a hash table ADT.
Definition htable.f90:36
Logging routines.
Definition log.f90:34
integer, parameter, public log_size
Definition log.f90:42
Definition math.f90:60
Defines a mesh.
Definition mesh.f90:34
subroutine mesh_generate_flags(this)
Definition mesh.f90:401
subroutine mesh_generate_facet_numbering(this)
Generate a unique facet numbering.
Definition mesh.f90:1133
subroutine mesh_init_dist(this, gdim, dist)
Initialise a mesh this based on a distribution dist.
Definition mesh.f90:205
logical function mesh_is_shared_point(this, p)
Check if a point is shared.
Definition mesh.f90:1833
integer, parameter, public neko_msh_max_zlbls
Max num. zone labels.
Definition mesh.f90:56
subroutine mesh_add_face(this, f)
Add a unique face represented as a 4-tuple to the mesh.
Definition mesh.f90:1473
integer function mesh_get_global_edge(this, e)
Return the global id of an edge e.
Definition mesh.f90:1783
subroutine mesh_generate_edge_conn(this)
Generate element-element connectivity via edges both between internal and between PEs.
Definition mesh.f90:889
subroutine mesh_add_edge(this, e)
Add a unique edge represented as a 2-tuple to the mesh.
Definition mesh.f90:1486
subroutine mesh_add_point(this, p, idx)
Add a unique point to the mesh.
Definition mesh.f90:1449
subroutine mesh_free(this)
Deallocate a mesh this.
Definition mesh.f90:312
subroutine mesh_mark_labeled_facet(this, f, e, label)
Mark facet f in element e with label.
Definition mesh.f90:1516
subroutine mesh_add_quad(this, el, el_glb, p1, p2, p3, p4)
Add a quadrilateral element to the mesh this.
Definition mesh.f90:1377
integer, parameter, public neko_msh_max_zlbl_len
Max length of a zone label.
Definition mesh.f90:58
subroutine mesh_mark_periodic_facet(this, f, e, pf, pe, pids)
Mark facet f in element e as periodic with (pf, pe)
Definition mesh.f90:1536
integer function mesh_get_local_edge(this, e)
Return the local id of an edge e.
Definition mesh.f90:1759
subroutine mesh_all_deformed(this)
Set all elements as if they are deformed.
Definition mesh.f90:434
logical function mesh_is_shared_edge(this, e)
Check if an edge is shared.
Definition mesh.f90:1847
subroutine mesh_init_common(this)
Definition mesh.f90:226
subroutine mesh_get_facet_ids(this, f, e, pids)
Get original ids of periodic points.
Definition mesh.f90:1550
subroutine mesh_create_periodic_ids(this, f, e, pf, pe)
Creates common ids for matching periodic points.
Definition mesh.f90:1626
subroutine mesh_reset_periodic_ids(this)
Reset ids of periodic points to their original ids.
Definition mesh.f90:1573
subroutine mesh_init_nelv(this, gdim, nelv)
Initialise a mesh this with nelv elements.
Definition mesh.f90:176
logical function mesh_is_shared_facet(this, f)
Check if a facet is shared.
Definition mesh.f90:1861
integer function mesh_get_global_facet(this, f)
Return the local id of a face f.
Definition mesh.f90:1803
integer function mesh_have_point_glb_idx(this, index)
Check if the mesh has a point given its global index.
Definition mesh.f90:1820
subroutine mesh_generate_external_point_conn(this)
Generate element-element connectivity via points between PEs.
Definition mesh.f90:821
subroutine mesh_apply_periodic_facet(this, f, e, pf, pe, pids)
Replaces the periodic point's id with a common id for matching periodic points.
Definition mesh.f90:1710
subroutine mesh_finalize(this)
Definition mesh.f90:386
subroutine mesh_generate_external_facet_conn(this)
Generate element-element connectivity via facets between PEs.
Definition mesh.f90:654
subroutine mesh_mark_curve_element(this, e, curve_data, curve_type)
Mark element e as a curve element.
Definition mesh.f90:1499
subroutine mesh_generate_conn(this)
Generate element-to-element connectivity.
Definition mesh.f90:440
subroutine mesh_add_hex(this, el, el_glb, p1, p2, p3, p4, p5, p6, p7, p8)
Add a hexahedral element to the mesh this.
Definition mesh.f90:1409
integer function mesh_get_local_facet(this, f)
Return the local id of a face f.
Definition mesh.f90:1771
integer function mesh_get_local_point(this, p)
Return the local id of a point p.
Definition mesh.f90:1742
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 a point.
Definition point.f90:35
Defines a quadrilateral element.
Definition quad.f90:34
integer, parameter, public neko_quad_neds
Number of edges.
Definition quad.f90:43
integer, parameter, public neko_quad_npts
Number of points.
Definition quad.f90:42
Implements a dynamic stack ADT.
Definition stack.f90:35
Implements a n-tuple.
Definition tuple.f90:34
Implements an unordered set ADT.
Definition uset.f90:35
Utilities.
Definition utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:266
Load-balanced linear distribution .
Definition datadist.f90:50
Base type for an element.
Definition element.f90:44
Hexahedron element.
Definition hex.f90:63
Integer based hash table.
Definition htable.f90:82
Integer 2-tuple based hash table.
Definition htable.f90:122
Integer 4-tuple based hash table.
Definition htable.f90:132
Integer*8 based hash table.
Definition htable.f90:92
Iterator for an integer based 2-tuple hash table.
Definition htable.f90:202
Iterator for an integer based 4-tuple hash table.
Definition htable.f90:211
Base type for a hash table.
Definition htable.f90:55
A point in with coordinates .
Definition point.f90:43
Quadrilateral element.
Definition quad.f90:58
Integer based stack.
Definition stack.f90:63
Integer 2-tuple based stack.
Definition stack.f90:84
Integer 4-tuple based stack.
Definition stack.f90:91
Integer*8 based stack.
Definition stack.f90:70
Integer based 4-tuple.
Definition tuple.f90:69
Integer based 2-tuple.
Definition tuple.f90:51
Integer*8 based unordered set.
Definition uset.f90:74
#define max(a, b)
Definition tensor.cu:40