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