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