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