Neko  0.9.99
A portable framework for high-order spectral element flow simulations
dofmap.f90
Go to the documentation of this file.
1 ! Copyright (c) 2020-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 !
35 module dofmap
36  use neko_config, only : neko_bcknd_device
37  use mesh, only : mesh_t
38  use space, only : space_t, gll
39  use tuple, only : tuple_i4_t, tuple4_i4_t
40  use num_types, only : i4, i8, rp, xp
41  use utils, only : neko_error, neko_warning
42  use fast3d, only : fd_weights_full
43  use tensor, only : tensr3, tnsr2d_el, trsp, addtnsr
44  use device
45  use math, only : add3, copy, rone, rzero
46  use element, only : element_t
47  use quad, only : quad_t
48  use hex, only : hex_t
49  use, intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
50  implicit none
51  private
52 
53  type, public :: dofmap_t
54  integer(kind=i8), allocatable :: dof(:,:,:,:)
55  logical, allocatable :: shared_dof(:,:,:,:)
56  real(kind=rp), allocatable :: x(:,:,:,:)
57  real(kind=rp), allocatable :: y(:,:,:,:)
58  real(kind=rp), allocatable :: z(:,:,:,:)
59  integer, private :: ntot
60 
61  type(mesh_t), pointer :: msh
62  type(space_t), pointer :: xh
63 
64  !
65  ! Device pointers (if present)
66  !
67  type(c_ptr) :: x_d = c_null_ptr
68  type(c_ptr) :: y_d = c_null_ptr
69  type(c_ptr) :: z_d = c_null_ptr
70 
71  contains
73  procedure, pass(this) :: init => dofmap_init
75  procedure, pass(this) :: free => dofmap_free
77  procedure, pass(this) :: size => dofmap_size
78  end type dofmap_t
79 
80 contains
81 
85  subroutine dofmap_init(this, msh, Xh)
86  class(dofmap_t) :: this
87  type(mesh_t), target, intent(inout) :: msh
88  type(space_t), target, intent(inout) :: Xh
89 
90  if ((msh%gdim .eq. 3 .and. xh%lz .eq. 1) .or. &
91  (msh%gdim .eq. 2 .and. xh%lz .gt. 1)) then
92  call neko_error("Invalid dimension of function space for the given mesh")
93  end if
94 
95  call this%free()
96 
97  this%msh => msh
98  this%Xh => xh
99 
100  this%ntot = xh%lx* xh%ly * xh%lz * msh%nelv
101 
102  !
103  ! Assign a unique id for all dofs
104  !
105 
106  allocate(this%dof(xh%lx, xh%ly, xh%lz, msh%nelv))
107  allocate(this%shared_dof(xh%lx, xh%ly, xh%lz, msh%nelv))
108 
109  this%dof = 0
110  this%shared_dof = .false.
111 
113  if (msh%gdim .eq. 3) then
114  call dofmap_number_points(this)
115  call dofmap_number_edges(this)
116  call dofmap_number_faces(this)
117  else
118  call dofmap_number_points(this)
119  call dofmap_number_edges(this)
120  end if
121 
122  !
123  ! Generate x,y,z-coordinates for all dofs
124  !
125 
126  allocate(this%x(xh%lx, xh%ly, xh%lz, msh%nelv))
127  allocate(this%y(xh%lx, xh%ly, xh%lz, msh%nelv))
128  allocate(this%z(xh%lx, xh%ly, xh%lz, msh%nelv))
129 
130  this%x = 0d0
131  this%y = 0d0
132  this%z = 0d0
134 
135  call dofmap_generate_xyz(this)
136 
137  if (neko_bcknd_device .eq. 1) then
138  call device_map(this%x, this%x_d, this%ntot)
139  call device_map(this%y, this%y_d, this%ntot)
140  call device_map(this%z, this%z_d, this%ntot)
141 
142  call device_memcpy(this%x, this%x_d, this%ntot, &
143  host_to_device, sync = .false.)
144  call device_memcpy(this%y, this%y_d, this%ntot, &
145  host_to_device, sync = .false.)
146  call device_memcpy(this%z, this%z_d, this%ntot, &
147  host_to_device, sync = .false.)
148  end if
149 
150  end subroutine dofmap_init
151 
153  subroutine dofmap_free(this)
154  class(dofmap_t), intent(inout) :: this
155 
156  if (allocated(this%dof)) then
157  deallocate(this%dof)
158  end if
159 
160  if (allocated(this%shared_dof)) then
161  deallocate(this%shared_dof)
162  end if
163 
164  if (allocated(this%x)) then
165  deallocate(this%x)
166  end if
167 
168  if (allocated(this%y)) then
169  deallocate(this%y)
170  end if
171 
172  if (allocated(this%z)) then
173  deallocate(this%z)
174  end if
175 
176  nullify(this%msh)
177  nullify(this%Xh)
178 
179  !
180  ! Cleanup the device (if present)
181  !
182  if (c_associated(this%x_d)) then
183  call device_free(this%x_d)
184  end if
185 
186  if (c_associated(this%y_d)) then
187  call device_free(this%y_d)
188  end if
189 
190  if (c_associated(this%z_d)) then
191  call device_free(this%z_d)
192  end if
193 
194  end subroutine dofmap_free
195 
197  pure function dofmap_size(this) result(res)
198  class(dofmap_t), intent(in) :: this
199  integer :: res
200  res = this%ntot
201  end function dofmap_size
202 
204  subroutine dofmap_number_points(this)
205  type(dofmap_t), target :: this
206  integer :: il, jl, ix, iy, iz
207  type(mesh_t), pointer :: msh
208  type(space_t), pointer :: Xh
209 
210  msh => this%msh
211  xh => this%Xh
212  do il = 1, msh%nelv
213  do jl = 1, msh%npts
214  ix = mod(jl - 1, 2) * (xh%lx - 1) + 1
215  iy = (mod(jl - 1, 4)/2) * (xh%ly - 1) + 1
216  iz = ((jl - 1)/4) * (xh%lz - 1) + 1
217  this%dof(ix, iy, iz, il) = int(msh%elements(il)%e%pts(jl)%p%id(), i8)
218  this%shared_dof(ix, iy, iz, il) = &
219  msh%is_shared(msh%elements(il)%e%pts(jl)%p)
220  end do
221  end do
222  end subroutine dofmap_number_points
223 
225  subroutine dofmap_number_edges(this)
226  type(dofmap_t), target :: this
227  type(mesh_t), pointer :: msh
228  type(space_t), pointer :: Xh
229  integer :: i,j,k
230  integer :: global_id
231  type(tuple_i4_t) :: edge
232  integer(kind=i8) :: num_dofs_edges(3) ! #dofs for each dir (r, s, t)
233  integer(kind=i8) :: edge_id, edge_offset
234  logical :: shared_dof
235 
236  msh => this%msh
237  xh => this%Xh
238 
239  ! Number of dofs on an edge excluding end-points
240  num_dofs_edges(1) = int(xh%lx - 2, i8)
241  num_dofs_edges(2) = int(xh%ly - 2, i8)
242  num_dofs_edges(3) = int(xh%lz - 2, i8)
243  edge_offset = int(msh%glb_mpts, i8) + int(1, 4)
244 
245  do i = 1, msh%nelv
246 
247  select type (ep => msh%elements(i)%e)
248  type is (hex_t)
249  !
250  ! Number edges in r-direction
251  !
252  call ep%edge_id(edge, 1)
253  shared_dof = msh%is_shared(edge)
254  global_id = msh%get_global(edge)
255  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
256  !Reverse order of tranversal if edge is reversed
257  if (int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) then
258  do concurrent(j = 2:xh%lx - 1)
259  k = xh%lx+1-j
260  this%dof(k, 1, 1, i) = edge_id + (j-2)
261  this%shared_dof(k, 1, 1, i) = shared_dof
262  end do
263  else
264  do concurrent(j = 2:xh%lx - 1)
265  k = j
266  this%dof(k, 1, 1, i) = edge_id + (j-2)
267  this%shared_dof(k, 1, 1, i) = shared_dof
268  end do
269  end if
270 
271  call ep%edge_id(edge, 3)
272  shared_dof = msh%is_shared(edge)
273  global_id = msh%get_global(edge)
274  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
275  if (int(edge%x(1), i8) .ne. this%dof(1, 1, xh%lz, i)) then
276  do concurrent(j = 2:xh%lx - 1)
277  k = xh%lx+1-j
278  this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
279  this%shared_dof(k, 1, xh%lz, i) = shared_dof
280  end do
281  else
282  do concurrent(j = 2:xh%lx - 1)
283  k = j
284  this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
285  this%shared_dof(k, 1, xh%lz, i) = shared_dof
286  end do
287  end if
288 
289  call ep%edge_id(edge, 2)
290  shared_dof = msh%is_shared(edge)
291  global_id = msh%get_global(edge)
292  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
293  if (int(edge%x(1), i8) .ne. this%dof(1, xh%ly, 1, i)) then
294  do concurrent(j = 2:xh%lx - 1)
295  k = xh%lx+1-j
296  this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
297  this%shared_dof(k, xh%ly, 1, i) = shared_dof
298  end do
299  else
300  do concurrent(j = 2:xh%lx - 1)
301  k = j
302  this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
303  this%shared_dof(k, xh%ly, 1, i) = shared_dof
304  end do
305  end if
306 
307  call ep%edge_id(edge, 4)
308  shared_dof = msh%is_shared(edge)
309  global_id = msh%get_global(edge)
310  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
311  if (int(edge%x(1), i8) .ne. this%dof(1, xh%ly, xh%lz, i)) then
312  do concurrent(j = 2:xh%lx - 1)
313  k = xh%lx+1-j
314  this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
315  this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
316  end do
317  else
318  do concurrent(j = 2:xh%lx - 1)
319  k = j
320  this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
321  this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
322  end do
323  end if
324 
325 
326  !
327  ! Number edges in s-direction
328  !
329  call ep%edge_id(edge, 5)
330  shared_dof = msh%is_shared(edge)
331  global_id = msh%get_global(edge)
332  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
333  if (int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) then
334  do concurrent(j = 2:xh%ly - 1)
335  k = xh%ly+1-j
336  this%dof(1, k, 1, i) = edge_id + (j-2)
337  this%shared_dof(1, k, 1, i) = shared_dof
338  end do
339  else
340  do concurrent(j = 2:xh%ly - 1)
341  k = j
342  this%dof(1, k, 1, i) = edge_id + (j-2)
343  this%shared_dof(1, k, 1, i) = shared_dof
344  end do
345  end if
346 
347  call ep%edge_id(edge, 7)
348  shared_dof = msh%is_shared(edge)
349  global_id = msh%get_global(edge)
350  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
351  if (int(edge%x(1), i8) .ne. this%dof(1, 1, xh%lz, i)) then
352  do concurrent(j = 2:xh%ly - 1)
353  k = xh%ly+1-j
354  this%dof(1, k, xh%lz, i) = edge_id + (j-2)
355  this%shared_dof(1, k, xh%lz, i) = shared_dof
356  end do
357  else
358  do concurrent(j = 2:xh%ly - 1)
359  k = j
360  this%dof(1, k, xh%lz, i) = edge_id + (j-2)
361  this%shared_dof(1, k, xh%lz, i) = shared_dof
362  end do
363  end if
364 
365  call ep%edge_id(edge, 6)
366  shared_dof = msh%is_shared(edge)
367  global_id = msh%get_global(edge)
368  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
369  if (int(edge%x(1), i8) .ne. this%dof(xh%lx, 1, 1, i)) then
370  do concurrent(j = 2:xh%ly - 1)
371  k = xh%ly+1-j
372  this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
373  this%shared_dof(xh%lx, k, 1, i) = shared_dof
374  end do
375  else
376  do concurrent(j = 2:xh%ly - 1)
377  k = j
378  this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
379  this%shared_dof(xh%lx, k, 1, i) = shared_dof
380  end do
381  end if
382 
383  call ep%edge_id(edge, i8)
384  shared_dof = msh%is_shared(edge)
385  global_id = msh%get_global(edge)
386  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
387  if (int(edge%x(1), i8) .ne. this%dof(xh%lx, 1, xh%lz, i)) then
388  do concurrent(j = 2:xh%ly - 1)
389  k = xh%lz+1-j
390  this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
391  this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
392  end do
393  else
394  do concurrent(j = 2:xh%ly - 1)
395  k = j
396  this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
397  this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
398  end do
399  end if
400 
401  !
402  ! Number edges in t-direction
403  !
404  call ep%edge_id(edge, 9)
405  shared_dof = msh%is_shared(edge)
406  global_id = msh%get_global(edge)
407  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
408  if (int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) then
409  do concurrent(j = 2:xh%lz - 1)
410  k = xh%lz+1-j
411  this%dof(1, 1, k, i) = edge_id + (j-2)
412  this%shared_dof(1, 1, k, i) = shared_dof
413  end do
414  else
415  do concurrent(j = 2:xh%lz - 1)
416  k = j
417  this%dof(1, 1, k, i) = edge_id + (j-2)
418  this%shared_dof(1, 1, k, i) = shared_dof
419  end do
420  end if
421 
422  call ep%edge_id(edge, 10)
423  shared_dof = msh%is_shared(edge)
424  global_id = msh%get_global(edge)
425  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
426  if (int(edge%x(1), i8) .ne. this%dof(xh%lx,1,1,i)) then
427  do concurrent(j = 2:xh%lz - 1)
428  k = xh%lz+1-j
429  this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
430  this%shared_dof(xh%lx, 1, k, i) = shared_dof
431  end do
432  else
433  do concurrent(j = 2:xh%lz - 1)
434  k = j
435  this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
436  this%shared_dof(xh%lx, 1, k, i) = shared_dof
437  end do
438  end if
439 
440  call ep%edge_id(edge, 11)
441  shared_dof = msh%is_shared(edge)
442  global_id = msh%get_global(edge)
443  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
444  if (int(edge%x(1), i8) .ne. this%dof(1, xh%ly, 1, i)) then
445  do concurrent(j = 2:xh%lz - 1)
446  k = xh%lz+1-j
447  this%dof(1, xh%ly, k, i) = edge_id + (j-2)
448  this%shared_dof(1, xh%ly, k, i) = shared_dof
449  end do
450  else
451  do concurrent(j = 2:xh%lz - 1)
452  k = j
453  this%dof(1, xh%ly, k, i) = edge_id + (j-2)
454  this%shared_dof(1, xh%ly, k, i) = shared_dof
455  end do
456  end if
457 
458  call ep%edge_id(edge, 12)
459  shared_dof = msh%is_shared(edge)
460  global_id = msh%get_global(edge)
461  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
462  if (int(edge%x(1), i8) .ne. this%dof(xh%lx, xh%ly, 1, i)) then
463  do concurrent(j = 2:xh%lz - 1)
464  k = xh%lz+1-j
465  this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
466  this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
467  end do
468  else
469  do concurrent(j = 2:xh%lz - 1)
470  k = j
471  this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
472  this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
473  end do
474  end if
475  type is (quad_t)
476  !
477  ! Number edges in r-direction
478  !
479  call ep%facet_id(edge, 3)
480  shared_dof = msh%is_shared(edge)
481  global_id = msh%get_global(edge)
482  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
483  !Reverse order of tranversal if edge is reversed
484  if (int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) then
485  do concurrent(j = 2:xh%lx - 1)
486  k = xh%lx+1-j
487  this%dof(k, 1, 1, i) = edge_id + (j-2)
488  this%shared_dof(k, 1, 1, i) = shared_dof
489  end do
490  else
491  do concurrent(j = 2:xh%lx - 1)
492  k = j
493  this%dof(k, 1, 1, i) = edge_id + (j-2)
494  this%shared_dof(k, 1, 1, i) = shared_dof
495  end do
496  end if
497 
498  call ep%facet_id(edge, 4)
499  shared_dof = msh%is_shared(edge)
500  global_id = msh%get_global(edge)
501  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
502  if (int(edge%x(1), i8) .ne. this%dof(1, xh%ly, 1, i)) then
503  do concurrent(j = 2:xh%lx - 1)
504  k = xh%lx+1-j
505  this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
506  this%shared_dof(k, xh%ly, 1, i) = shared_dof
507  end do
508  else
509  do concurrent(j = 2:xh%lx - 1)
510  k = j
511  this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
512  this%shared_dof(k, xh%ly, 1, i) = shared_dof
513  end do
514  end if
515 
516  !
517  ! Number edges in s-direction
518  !
519  call ep%facet_id(edge, 1)
520  shared_dof = msh%is_shared(edge)
521  global_id = msh%get_global(edge)
522  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
523  if (int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) then
524  do concurrent(j = 2:xh%ly - 1)
525  k = xh%ly+1-j
526  this%dof(1, k, 1, i) = edge_id + (j-2)
527  this%shared_dof(1, k, 1, i) = shared_dof
528  end do
529  else
530  do concurrent(j = 2:xh%ly - 1)
531  k = j
532  this%dof(1, k, 1, i) = edge_id + (j-2)
533  this%shared_dof(1, k, 1, i) = shared_dof
534  end do
535  end if
536 
537  call ep%facet_id(edge, 2)
538  shared_dof = msh%is_shared(edge)
539  global_id = msh%get_global(edge)
540  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
541  if (int(edge%x(1), i8) .ne. this%dof(xh%lx,1,1,i)) then
542  do concurrent(j = 2:xh%ly - 1)
543  k = xh%ly+1-j
544  this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
545  this%shared_dof(xh%lx, k, 1, i) = shared_dof
546  end do
547  else
548  do concurrent(j = 2:xh%ly - 1)
549  k = j
550  this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
551  this%shared_dof(xh%lx, k, 1, i) = shared_dof
552  end do
553  end if
554  end select
555 
556  end do
557  end subroutine dofmap_number_edges
558 
560  subroutine dofmap_number_faces(this)
561  type(dofmap_t), target :: this
562  type(mesh_t), pointer :: msh
563  type(space_t), pointer :: Xh
564  integer :: i,j,k
565  integer :: global_id
566  type(tuple4_i4_t) :: face, face_order
567  integer(kind=i8) :: num_dofs_faces(3) ! #dofs for each dir (r, s, t)
568  integer(kind=i8) :: facet_offset, facet_id
569  logical :: shared_dof
570 
571  msh => this%msh
572  xh => this%Xh
573 
575  facet_offset = int(msh%glb_mpts, i8) + &
576  int(msh%glb_meds, i8) * int(xh%lx-2, i8) + int(1, i8)
577 
578  ! Number of dofs on an face excluding end-points
579  num_dofs_faces(1) = int((xh%ly - 2) * (xh%lz - 2), i8)
580  num_dofs_faces(2) = int((xh%lx - 2) * (xh%lz - 2), i8)
581  num_dofs_faces(3) = int((xh%lx - 2) * (xh%ly - 2), i8)
582 
583  do i = 1, msh%nelv
584 
585  !
586  ! Number facets in r-direction (s, t)-plane
587  !
588  call msh%elements(i)%e%facet_id(face, 1)
589  call msh%elements(i)%e%facet_order(face_order, 1)
590  shared_dof = msh%is_shared(face)
591  global_id = msh%get_global(face)
592  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(1)
593  do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
594  this%dof(1, j, k, i) = &
595  dofmap_facetidx(face_order, face, facet_id, j, k, xh%lz, xh%ly)
596  this%shared_dof(1, j, k, i) = shared_dof
597  end do
598 
599  call msh%elements(i)%e%facet_id(face, 2)
600  call msh%elements(i)%e%facet_order(face_order, 2)
601  shared_dof = msh%is_shared(face)
602  global_id = msh%get_global(face)
603  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(1)
604  do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
605  this%dof(xh%lx, j, k, i) = &
606  dofmap_facetidx(face_order, face, facet_id, j, k, xh%lz, xh%ly)
607  this%shared_dof(xh%lx, j, k, i) = shared_dof
608  end do
609 
610 
611  !
612  ! Number facets in s-direction (r, t)-plane
613  !
614  call msh%elements(i)%e%facet_id(face, 3)
615  call msh%elements(i)%e%facet_order(face_order, 3)
616  shared_dof = msh%is_shared(face)
617  global_id = msh%get_global(face)
618  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(2)
619  do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
620  this%dof(j, 1, k, i) = &
621  dofmap_facetidx(face_order, face, facet_id, k, j, xh%lz, xh%lx)
622  this%shared_dof(j, 1, k, i) = shared_dof
623  end do
624 
625  call msh%elements(i)%e%facet_id(face, 4)
626  call msh%elements(i)%e%facet_order(face_order, 4)
627  shared_dof = msh%is_shared(face)
628  global_id = msh%get_global(face)
629  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(2)
630  do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
631  this%dof(j, xh%ly, k, i) = &
632  dofmap_facetidx(face_order, face, facet_id, k, j, xh%lz, xh%lx)
633  this%shared_dof(j, xh%ly, k, i) = shared_dof
634  end do
635 
636 
637  !
638  ! Number facets in t-direction (r, s)-plane
639  !
640  call msh%elements(i)%e%facet_id(face, 5)
641  call msh%elements(i)%e%facet_order(face_order, 5)
642  shared_dof = msh%is_shared(face)
643  global_id = msh%get_global(face)
644  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(3)
645  do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
646  this%dof(j, k, 1, i) = &
647  dofmap_facetidx(face_order, face, facet_id, k, j, xh%ly, xh%lx)
648  this%shared_dof(j, k, 1, i) = shared_dof
649  end do
650 
651  call msh%elements(i)%e%facet_id(face, 6)
652  call msh%elements(i)%e%facet_order(face_order, 6)
653  shared_dof = msh%is_shared(face)
654  global_id = msh%get_global(face)
655  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(3)
656  do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
657  this%dof(j, k, xh%lz, i) = &
658  dofmap_facetidx(face_order, face, facet_id, k, j, xh%lz, xh%lx)
659  this%shared_dof(j, k, xh%lz, i) = shared_dof
660  end do
661  end do
662 
663  end subroutine dofmap_number_faces
664 
666  pure function dofmap_facetidx(face_order, face, facet_id, k1, j1, lk1, &
667  lj1) result(facet_idx)
668  type(tuple4_i4_t), intent(in) :: face_order, face
669  integer(kind=i8), intent(in) :: facet_id
670  integer(kind=i8) :: facet_idx
671  integer, intent(in) :: k1, j1, lk1, lj1
672  integer :: k, j, lk, lj
673 
674  k = k1 - 2
675  j = j1 - 2
676  lk = lk1 - 2
677  lj = lj1 - 2
678 
679  ! Given the indexes k,j for a GLL point on the inner part of the
680  ! face, we assign a unique number to it that depends on the
681  ! corner with the lowest id and its neighbour with the lowest
682  ! id. The id is assigned in this way to be consistent regardless
683  ! of how the faces are rotated or mirrored.
684  !
685  ! 4 -------- 3
686  ! | | k
687  ! |----->| ^
688  ! |----->| |
689  ! |----->| |
690  ! 1 -------- 2 0--->j
691 
692 
693  if (face_order%x(1) .eq. face%x(1)) then
694  if (face_order%x(2) .lt. face_order%x(4)) then
695  facet_idx = facet_id + j + k*lj
696  else
697  facet_idx = facet_id + j*lk + k
698  end if
699  else if (face_order%x(2) .eq. face%x(1)) then
700  if (face_order%x(3) .lt. face_order%x(1)) then
701  facet_idx = facet_id + lk*(lj-1-j) + k
702  else
703  facet_idx = facet_id + (lj-1-j) + k*lj
704  end if
705  else if (face_order%x(3) .eq. face%x(1)) then
706  if (face_order%x(4) .lt. face_order%x(2)) then
707  facet_idx = facet_id + (lj-1-j) + lj*(lk-1-k)
708  else
709  facet_idx = facet_id + lk*(lj-1-j) + (lk-1-k)
710  end if
711  else if (face_order%x(4) .eq. face%x(1)) then
712  if (face_order%x(1) .lt. face_order%x(3)) then
713  facet_idx = facet_id + lk*j + (lk-1-k)
714  else
715  facet_idx = facet_id + j + lj*(lk-1-k)
716  end if
717  end if
718 
719  end function dofmap_facetidx
720 
723  subroutine dofmap_generate_xyz(this)
724  type(dofmap_t), target :: this
725  integer :: i, j, el_idx
726  type(mesh_t), pointer :: msh
727  type(space_t), pointer :: Xh
728  real(kind=rp) :: rp_curve_data(5), curve_data_tot(5,12)
729  logical :: midpoint
730  integer :: n_edge, curve_type(12)
731 
732  msh => this%msh
733  xh => this%Xh
734 
735  if (msh%gdim .eq. 3) then
736  n_edge = 12
737  else
738  n_edge = 4
739  end if
740 
741  do i = 1, msh%nelv
742  call dofmap_xyzlin(xh, msh, msh%elements(i)%e, this%x(1,1,1,i), &
743  this%y(1,1,1,i), this%z(1,1,1,i))
744  end do
745  do i = 1, msh%curve%size
746  midpoint = .false.
747  el_idx = msh%curve%curve_el(i)%el_idx
748  curve_type = msh%curve%curve_el(i)%curve_type
749  curve_data_tot = msh%curve%curve_el(i)%curve_data
750  do j = 1, n_edge
751  if (curve_type(j) .eq. 4) then
752  midpoint = .true.
753  end if
754  end do
755  if (midpoint .and. xh%lx .gt. 2) then
756  call dofmap_xyzquad(xh, msh, msh%elements(el_idx)%e, &
757  this%x(1, 1, 1, el_idx), this%y(1, 1, 1, el_idx), &
758  this%z(1 ,1, 1, el_idx), curve_type, curve_data_tot)
759  end if
760  end do
761  do i = 1, msh%curve%size
762  el_idx = msh%curve%curve_el(i)%el_idx
763  do j = 1, 8
764  if (msh%curve%curve_el(i)%curve_type(j) .eq. 3) then
765  rp_curve_data = msh%curve%curve_el(i)%curve_data(1:5,j)
766  call arc_surface(j, rp_curve_data, &
767  this%x(1, 1, 1, el_idx), &
768  this%y(1, 1, 1, el_idx), &
769  this%z(1, 1, 1, el_idx), &
770  xh, msh%elements(el_idx)%e, msh%gdim)
771  end if
772  end do
773  end do
774  if (associated(msh%apply_deform)) then
775  call msh%apply_deform(this%x, this%y, this%z, xh%lx, xh%ly, xh%lz)
776  end if
777  end subroutine dofmap_generate_xyz
778 
787  subroutine dofmap_xyzlin(Xh, msh, element, x, y, z)
788  type(mesh_t), pointer, intent(in) :: msh
789  type(space_t), intent(in) :: Xh
790  class(element_t), intent(in) :: element
791  real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz), &
792  y(xh%lx, xh%ly, xh%lz), &
793  z(xh%lx, xh%ly, xh%lz)
794  real(kind=rp) :: xyzb(2,2,2,3), zgml(xh%lx, 3)
795  real(kind=rp) :: jx(xh%lx*2)
796  real(kind=rp) :: jxt(xh%lx*2), jyt(xh%lx*2), jzt(xh%lx*2)
797  real(kind=rp) :: w(4*xh%lx**3), tmp(xh%lx, xh%lx, xh%lx)
798  real(kind=rp), dimension(2), parameter :: zlin = [-1d0, 1d0]
799 
800  integer :: j, k
801 
802  zgml = 0d0
803  xyzb = 0d0
804 
805  w = 0d0
806  call copy(zgml(1,1), xh%zg(1,1), xh%lx)
807  call copy(zgml(1,2), xh%zg(1,2), xh%ly)
808  if (msh%gdim .gt. 2) then
809  call copy(zgml(1,3), xh%zg(1,3), xh%lz)
810  end if
811 
812  k = 1
813  do j = 1, xh%lx
814  call fd_weights_full(zgml(j,1), zlin, 1, 0, jxt(k))
815  call fd_weights_full(zgml(j,2), zlin, 1, 0, jyt(k))
816  if (msh%gdim .gt. 2) then
817  call fd_weights_full(zgml(j,3), zlin, 1, 0, jzt(k))
818  end if
819  k = k + 2
820  end do
821  call trsp(jx, xh%lx, jxt, 2)
822 
823  if (msh%gdim .eq. 2) then
824  jzt = 1d0
825  end if
826 
827  if (msh%gdim .gt. 2) then
828  do concurrent(j = 1:msh%gdim)
829  xyzb(1,1,1,j) = element%pts(1)%p%x(j)
830  xyzb(2,1,1,j) = element%pts(2)%p%x(j)
831  xyzb(1,2,1,j) = element%pts(3)%p%x(j)
832  xyzb(2,2,1,j) = element%pts(4)%p%x(j)
833 
834  xyzb(1,1,2,j) = element%pts(5)%p%x(j)
835  xyzb(2,1,2,j) = element%pts(6)%p%x(j)
836  xyzb(1,2,2,j) = element%pts(7)%p%x(j)
837  xyzb(2,2,2,j) = element%pts(8)%p%x(j)
838  end do
839  else
840  do concurrent(j = 1:msh%gdim)
841  xyzb(1,1,1,j) = element%pts(1)%p%x(j)
842  xyzb(2,1,1,j) = element%pts(2)%p%x(j)
843  xyzb(1,2,1,j) = element%pts(3)%p%x(j)
844  xyzb(2,2,1,j) = element%pts(4)%p%x(j)
845  end do
846  end if
847  if (msh%gdim .eq. 3) then
848  call tensr3(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt, jzt, w)
849  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
850  call tensr3(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt, jzt, w)
851  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
852  call tensr3(tmp, xh%lz, xyzb(1,1,1,3), 2, jx, jyt, jzt, w)
853  call copy(z, tmp, xh%lx*xh%ly*xh%lz)
854  else
855  call tnsr2d_el(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt)
856  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
857  call tnsr2d_el(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt)
858  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
859  end if
860  end subroutine dofmap_xyzlin
861 
862  subroutine dofmap_xyzquad(Xh, msh, element, x, y, z, curve_type, curve_data)
863  type(mesh_t), pointer, intent(in) :: msh
864  type(space_t), intent(in) :: Xh
865  class(element_t), intent(in) :: element
866  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz), intent(inout) :: x, y, z
867  integer :: curve_type(12), eindx(12)
868  real(kind=rp) :: curve_data(5,12), x3(3,3,3), y3(3,3,3), z3(3,3,3)
869  type(space_t), target :: xh3
870  real(kind=rp), dimension(3), parameter :: zquad = [-1d0, 0d0,1d0]
871  real(kind=rp) :: zg(3)
872  real(kind=rp), dimension(Xh%lx, Xh%lx, Xh%lx) :: tmp
873  real(kind=rp) :: jx(xh%lx*3)
874  real(kind=rp) :: jxt(xh%lx*3), jyt(xh%lx*3), jzt(xh%lx*3)
875  real(kind=rp) :: w(4*xh%lxyz,2)
876  integer :: j, k, n_edges
877  eindx = [2 , 6 , 8 , 4, &
878  20 , 24 , 26 , 22, &
879  10 , 12 , 18 , 16]
880 
881  w = 0d0
882  if (msh%gdim .eq. 3) then
883  n_edges = 12
884  call xh3%init(gll, 3, 3, 3)
885  else
886  n_edges = 4
887  call xh3%init(gll, 3, 3)
888  end if
889  call dofmap_xyzlin(xh3, msh, element, x3, y3, z3)
890 
891  do k = 1, n_edges
892  if (curve_type(k) .eq. 4) then
893  x3(eindx(k),1,1) = curve_data(1,k)
894  y3(eindx(k),1,1) = curve_data(2,k)
895  z3(eindx(k),1,1) = curve_data(3,k)
896  end if
897  end do
898  zg(1) = -1
899  zg(2) = 0
900  zg(3) = 1
901  if (msh%gdim .eq. 3) then
902  call gh_face_extend_3d(x3, zg, 3, 2, w(1,1), w(1,2)) ! 2 --> edge extend
903  call gh_face_extend_3d(y3, zg, 3, 2, w(1,1), w(1,2))
904  call gh_face_extend_3d(z3, zg, 3, 2, w(1,1), w(1,2))
905  else
906  call neko_warning(' m deformation not supported for 2d yet')
907  call gh_face_extend_2d(x3, zg, 3, 2, w(1,1), w(1,2)) ! 2 --> edge extend
908  call gh_face_extend_2d(y3, zg, 3, 2, w(1,1), w(1,2))
909  end if
910  k = 1
911  do j = 1, xh%lx
912  call fd_weights_full(xh%zg(j,1), zquad, 2, 0, jxt(k))
913  call fd_weights_full(xh%zg(j,2), zquad, 2, 0, jyt(k))
914  if (msh%gdim .gt. 2) then
915  call fd_weights_full(xh%zg(j,3), zquad, 2, 0, jzt(k))
916  end if
917  k = k + 3
918  end do
919  call trsp(jx, xh%lx, jxt, 3)
920  if (msh%gdim .eq. 3) then
921  call tensr3(tmp, xh%lx, x3, 3, jx, jyt, jzt, w)
922  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
923  call tensr3(tmp, xh%ly, y3, 3, jx, jyt, jzt, w)
924  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
925  call tensr3(tmp, xh%lz, z3, 3, jx, jyt, jzt, w)
926  call copy(z, tmp, xh%lx*xh%ly*xh%lz)
927  else
928  call tnsr2d_el(tmp, xh%lx, x3, 3, jx, jyt)
929  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
930  call tnsr2d_el(tmp, xh%ly, y3, 3, jx, jyt)
931  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
932  end if
933 
934  call xh3%free()
935  end subroutine dofmap_xyzquad
936 
942  subroutine gh_face_extend_3d(x, zg, n, gh_type, e, v)
943  integer, intent(in) :: n
944  real(kind=rp), intent(inout) :: x(n, n, n)
945  real(kind=rp), intent(in) :: zg(n)
946  real(kind=rp), intent(inout) :: e(n, n, n)
947  real(kind=rp), intent(inout) :: v(n, n, n)
948  integer :: gh_type, ntot, kk, jj, ii, k, j, i
949  real(kind=xp) :: si, sj, sk, hi, hj, hk
950 
951  !
952  ! Build vertex interpolant
953  !
954  ntot = n**3
955  do concurrent(i = 1:ntot)
956  v(i,1,1) = 0.0_rp
957  end do
958 
959  do concurrent(i = 1:n, j = 1:n, k = 1:n, &
960  ii = 1:n:n-1, jj = 1:n:n-1, kk = 1:n:n-1)
961  si = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
962  sj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
963  sk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
964  v(i,j,k) = v(i,j,k) + si * sj* sk * x(ii, jj, kk)
965  end do
966 
967  if (gh_type .eq. 1) then
968  do concurrent(i = 1:ntot)
969  x(i,1,1) = v(i,1,1)
970  end do
971  return
972  end if
973  !
974  !
975  ! Extend 12 edges
976  do concurrent(i = 1:ntot)
977  e(i,1,1) = 0.0_rp
978  end do
979  !
980  ! x-edges
981  !
982  do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1, kk = 1:n:n-1)
983  hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
984  hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
985  e(i,j,k) = e(i,j,k) + hj*hk*(x(i, jj, kk) - v(i, jj, kk))
986  end do
987  !
988  ! y-edges
989  !
990  do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, kk = 1:n:n-1)
991  hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
992  hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
993  e(i,j,k) = e(i,j,k) + hi*hk*(x(ii, j, kk) - v(ii, j, kk))
994  end do
995  !
996  ! z-edges
997  !
998  do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, jj = 1:n:n-1)
999  hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1000  hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1001  e(i,j,k) = e(i,j,k) + hi*hj*(x(ii, jj, k) - v(ii, jj, k))
1002  end do
1003 
1004  do concurrent(i = 1:ntot)
1005  e(i,1,1) = e(i,1,1) + v(i,1,1)
1006  end do
1007 
1008  if (gh_type .eq. 2) then
1009  do concurrent(i = 1:ntot)
1010  x(i,1,1) = e(i,1,1)
1011  end do
1012  return
1013  end if
1014  !
1015  ! Extend faces
1016  !
1017  do concurrent(i = 1:ntot)
1018  v(i,1,1) = 0.0_rp
1019  end do
1020  !
1021  ! x-edges
1022  !
1023  do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1)
1024  hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1025  v(i,j,k) = v(i,j,k) + hi*(x(ii,j,k)-e(ii,j,k))
1026  end do
1027 
1028  !
1029  ! y-edges
1030  !
1031  do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1)
1032  hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1033  v(i,j,k) = v(i,j,k) + hj*(x(i, jj, k) - e(i, jj, k))
1034  end do
1035 
1036  !
1037  ! z-edges
1038  !
1039  do concurrent(i = 1:n, j = 1:n, k = 1:n, kk = 1:n:n-1)
1040  hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
1041  v(i,j,k) = v(i,j,k) + hk*(x(i, j, kk) - e(i, j, kk))
1042  end do
1043 
1044  do concurrent(i = 1:ntot)
1045  v(i,1,1) = v(i,1,1) + e(i,1,1)
1046  x(i,1,1) = v(i,1,1)
1047  end do
1048 
1049  end subroutine gh_face_extend_3d
1050 
1054  subroutine gh_face_extend_2d(x, zg, n, gh_type, e, v)
1055  integer, intent(in) :: n
1056  real(kind=rp), intent(inout) :: x(n, n)
1057  real(kind=rp), intent(in) :: zg(n)
1058  real(kind=rp), intent(inout) :: e(n, n)
1059  real(kind=rp), intent(inout) :: v(n, n)
1060  integer, intent(in) :: gh_type
1061  integer :: i,j , jj, ii, ntot
1062  real(kind=rp) :: si, sj, hi, hj
1063 
1064  !Build vertex interpolant
1065 
1066  ntot = n*n
1067  call rzero(v, ntot)
1068  do jj = 1, n, n-1
1069  do ii = 1, n, n-1
1070  do j = 1, n
1071  do i = 1, n
1072  si = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1073  sj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1074  v(i,j) = v(i,j) + si*sj*x(ii, jj)
1075  end do
1076  end do
1077  end do
1078  end do
1079  if (gh_type .eq. 1) then
1080  call copy(x, v, ntot)
1081  return
1082  end if
1083 
1084  !Extend 4 edges
1085  call rzero(e, ntot)
1086 
1087  !x-edges
1088 
1089  do jj = 1, n, n-1
1090  do j = 1, n
1091  do i = 1, n
1092  hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1093  e(i,j) = e(i,j) + hj*(x(i, jj) - v(i, jj))
1094  end do
1095  end do
1096  end do
1097 
1098  !y-edges
1099 
1100  do ii = 1, n, n-1
1101  do j = 1, n
1102  do i = 1, n
1103  hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1104  e(i,j) = e(i,j) + hi*(x(ii,j)-v(ii,j))
1105  end do
1106  end do
1107  end do
1108 
1109  call add3(x, e, v, ntot)
1110 
1111  end subroutine gh_face_extend_2d
1112 
1113 
1114 
1115  subroutine arc_surface(isid, curve_data, x, y, z, Xh, element, gdim)
1116  integer, intent(in) :: isid, gdim
1117  type(space_t), intent(in) :: Xh
1118  class(element_t) :: element
1119  real(kind=rp), dimension(5), intent(in) :: curve_data
1120  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz), intent(inout) :: x, y, z
1121  real(kind=rp) :: pt1x, pt1y, pt2x, pt2y, pt12x, pt12y
1122  real(kind=rp) :: radius, dtheta, r, xys
1123  real(kind=rp) :: theta0, xcenn, ycenn, h(xh%lx, 3, 2)
1124  real(kind=rp) :: xcrved(xh%lx), ycrved(xh%lx), xs, ys
1125  integer :: isid1, ixt, iyt, izt, ix, itmp
1126  ! Cyclic to symmetric face mapping
1127  integer(i4), dimension(6), parameter :: fcyc_to_sym = [3, 2, 4, 1, 5, 6]
1128  ! Cyclic to symmetric edge mapping
1129  integer(i4), dimension(12), parameter :: ecyc_to_sym = [1, 6, 2, 5, 3, 8,&
1130  & 4, 7, 9, 10, 12, 11]
1131  ! Symmetric edge to vertex mapping
1132  integer, parameter, dimension(2, 12) :: edge_nodes = reshape([1, 2, 3, 4, &
1133  & 5, 6, 7, 8, 1, 3, 2, 4, 5, 7, 6, 8, 1, 5, 2, 6, 3, 7, 4, 8], &
1134  & [2,12])
1135  ! copy from hex as this has private attribute there
1136 
1137  ! this subroutine is a mess of symmetric and cyclic edge/face numberring and
1138  ! cannot be cleaned without changing an input format (isid seems to be
1139  ! a cyclic edge number)
1140  ! following according to cyclic edge numbering and orientation
1141  itmp = ecyc_to_sym(isid)
1142  select case (isid)
1143  case (1:2,5:6)
1144  pt1x = element%pts(edge_nodes(1, itmp))%p%x(1)
1145  pt1y = element%pts(edge_nodes(1, itmp))%p%x(2)
1146  pt2x = element%pts(edge_nodes(2, itmp))%p%x(1)
1147  pt2y = element%pts(edge_nodes(2, itmp))%p%x(2)
1148  case (3:4,7:8)
1149  pt1x = element%pts(edge_nodes(2, itmp))%p%x(1)
1150  pt1y = element%pts(edge_nodes(2, itmp))%p%x(2)
1151  pt2x = element%pts(edge_nodes(1, itmp))%p%x(1)
1152  pt2y = element%pts(edge_nodes(1, itmp))%p%x(2)
1153  end select
1154  ! find slope of perpendicular
1155  radius = curve_data(1)
1156  xs = pt2y-pt1y
1157  ys = pt1x-pt2x
1158  ! make length radius
1159  xys = sqrt(xs**2 + ys**2)
1160  ! sanity check
1161  if (abs(2.0 * radius) <= xys * 1.00001) &
1162  & call neko_error('Radius to small for arced element surface')
1163  ! find center
1164  dtheta = abs(asin(0.5_xp*xys/radius))
1165  pt12x = (pt1x + pt2x)/2.0
1166  pt12y = (pt1y + pt2y)/2.0
1167  xcenn = pt12x - xs/xys * radius*cos(dtheta)
1168  ycenn = pt12y - ys/xys * radius*cos(dtheta)
1169  theta0 = atan2((pt12y-ycenn), (pt12x-xcenn))
1170 ! compute perturbation of geometry
1171  isid1 = mod(isid+4-1, 4)+1
1172  call compute_h(h, xh%zg, gdim, xh%lx)
1173  if (radius < 0.0) dtheta = -dtheta
1174  do ix = 1, xh%lx
1175  ixt = ix
1176  if (isid1 .gt. 2) ixt = xh%lx+1-ix
1177  r = xh%zg(ix,1)
1178  xcrved(ixt) = xcenn + abs(radius) * cos(theta0 + r*dtheta) &
1179  - ( h(ix,1,1)*pt1x + h(ix,1,2)*pt2x )
1180  ycrved(ixt) = ycenn + abs(radius) * sin(theta0 + r*dtheta) &
1181  - ( h(ix,1,1)*pt1y + h(ix,1,2)*pt2y )
1182  end do
1183 ! points all set, add perturbation to current mesh.
1184 ! LEGACY WARNING
1185 ! I dont want to dive in this again, Martin Karp 2/3 - 2021
1186  isid1 = fcyc_to_sym(isid1)
1187  izt = (isid-1)/4+1
1188  iyt = isid1-2
1189  ixt = isid1
1190  if (isid1 .le. 2) then
1191  call addtnsr(x, h(1, 1, ixt), xcrved, h(1, 3, izt), &
1192  xh%lx, xh%ly, xh%lz)
1193  call addtnsr(y, h(1, 1, ixt), ycrved, h(1, 3, izt), &
1194  xh%lx, xh%ly, xh%lz)
1195  else
1196  call addtnsr(x, xcrved, h(1, 2, iyt), h(1, 3, izt), &
1197  xh%lx, xh%ly, xh%lz)
1198  call addtnsr(y, ycrved, h(1, 2, iyt), h(1, 3, izt), &
1199  xh%lx, xh%ly, xh%lz)
1200  end if
1201  end subroutine arc_surface
1202 
1203  subroutine compute_h(h, zgml, gdim, lx)
1204  integer, intent(in) :: lx, gdim
1205  real(kind=rp), intent(inout) :: h(lx, 3, 2)
1206  real(kind=rp), intent(in) :: zgml(lx, 3)
1207  integer :: ix, iy, iz
1208 
1209  do ix = 1, lx
1210  h(ix,1,1) = (1.0_rp - zgml(ix, 1)) * 0.5_rp
1211  h(ix,1,2) = (1.0_rp + zgml(ix, 1)) * 0.5_rp
1212  end do
1213 
1214  do iy = 1, lx
1215  h(iy,2,1) = (1.0_rp - zgml(iy, 2)) * 0.5_rp
1216  h(iy,2,2) = (1.0_rp + zgml(iy, 2)) * 0.5_rp
1217  end do
1218 
1219  if (gdim .eq. 3) then
1220  do iz = 1, lx
1221  h(iz,3,1) = (1.0_rp - zgml(iz, 3)) * 0.5_rp
1222  h(iz,3,2) = (1.0_rp + zgml(iz, 3)) * 0.5_rp
1223  end do
1224  else
1225  call rone(h(1,3,1), lx)
1226  call rone(h(1,3,2), lx)
1227  end if
1228 
1229  end subroutine compute_h
1230 
1231 end module dofmap
Map a Fortran array to a device (allocate and associate)
Definition: device.F90:57
Copy data between host and device (or device and device)
Definition: device.F90:51
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
subroutine, public device_free(x_d)
Deallocate memory on the device.
Definition: device.F90:185
Defines a mapping of the degrees of freedom.
Definition: dofmap.f90:35
subroutine gh_face_extend_2d(x, zg, n, gh_type, e, v)
Extend 2D faces into interior via gordon hall gh_type: 1 - vertex only 2 - vertex and faces.
Definition: dofmap.f90:1055
subroutine dofmap_xyzquad(Xh, msh, element, x, y, z, curve_type, curve_data)
Definition: dofmap.f90:863
subroutine dofmap_generate_xyz(this)
Generate x,y,z-coordinates for all dofs.
Definition: dofmap.f90:724
subroutine compute_h(h, zgml, gdim, lx)
Definition: dofmap.f90:1204
subroutine arc_surface(isid, curve_data, x, y, z, Xh, element, gdim)
Definition: dofmap.f90:1116
subroutine dofmap_free(this)
Destructor.
Definition: dofmap.f90:154
subroutine dofmap_init(this, msh, Xh)
Constructor.
Definition: dofmap.f90:86
subroutine dofmap_number_edges(this)
Assing numbers to dofs on edges.
Definition: dofmap.f90:226
subroutine dofmap_number_points(this)
Assign numbers to each dofs on points.
Definition: dofmap.f90:205
subroutine dofmap_xyzlin(Xh, msh, element, x, y, z)
Generate the x, y, z coordinates of the dofs in a signle element, assuming linear element edges.
Definition: dofmap.f90:788
pure integer function dofmap_size(this)
Return the total number of dofs in the dofmap, lx*ly*lz*nelv.
Definition: dofmap.f90:198
subroutine dofmap_number_faces(this)
Assign numbers to dofs on faces.
Definition: dofmap.f90:561
subroutine gh_face_extend_3d(x, zg, n, gh_type, e, v)
Extend faces into interior via gordon hall gh_type: 1 - vertex only 2 - vertex and edges 3 - vertex,...
Definition: dofmap.f90:943
pure integer(kind=i8) function dofmap_facetidx(face_order, face, facet_id, k1, j1, lk1, lj1)
Get idx for GLL point on face depending on face ordering k and j.
Definition: dofmap.f90:668
Fast diagonalization methods from NEKTON.
Definition: fast3d.f90:61
subroutine, public fd_weights_full(xi, x, n, m, c)
Compute finite-difference stencil weights for evaluating derivatives up to order at a point.
Definition: fast3d.f90:105
Defines a hexahedron element.
Definition: hex.f90:34
Definition: math.f90:60
subroutine, public rone(a, n)
Set all elements to one.
Definition: math.f90:228
subroutine, public add3(a, b, c, n)
Vector addition .
Definition: math.f90:600
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:239
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:195
Defines a mesh.
Definition: mesh.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
integer, parameter, public i8
Definition: num_types.f90:7
integer, parameter, public i4
Definition: num_types.f90:6
integer, parameter, public xp
Definition: num_types.f90:14
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines a quadrilateral element.
Definition: quad.f90:34
Defines a function space.
Definition: space.f90:34
integer, parameter, public gll
Definition: space.f90:48
Tensor operations.
Definition: tensor.f90:61
subroutine, public tensr3(v, nv, u, nu, A, Bt, Ct, w)
Tensor product .
Definition: tensor.f90:91
subroutine, public addtnsr(s, h1, h2, h3, nx, ny, nz)
Maps and adds to S a tensor product form of the three functions H1,H2,H3. This is a single element ro...
Definition: tensor.f90:265
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
Definition: tensor.f90:121
subroutine, public tnsr2d_el(v, nv, u, nu, A, Bt)
Computes .
Definition: tensor.f90:154
Implements a n-tuple.
Definition: tuple.f90:34
Utilities.
Definition: utils.f90:35
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition: utils.f90:245
Base type for an element.
Definition: element.f90:44
Hexahedron element.
Definition: hex.f90:63
Quadrilateral element.
Definition: quad.f90:58
The function space for the SEM solution fields.
Definition: space.f90:62
Integer based 4-tuple.
Definition: tuple.f90:69
Integer based 2-tuple.
Definition: tuple.f90:51