Neko  0.8.1
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  do j = 2, xh%lx - 1
256  k = j
257  if(int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) k = xh%lx+1-j
258  this%dof(k, 1, 1, i) = edge_id
259  this%shared_dof(k, 1, 1, i) = shared_dof
260  edge_id = edge_id + 1
261  end do
262 
263  call ep%edge_id(edge, 3)
264  shared_dof = msh%is_shared(edge)
265  global_id = msh%get_global(edge)
266  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
267  do j = 2, xh%lx - 1
268  k = j
269  if(int(edge%x(1), i8) .ne. this%dof(1,1,xh%lz,i)) k = xh%lx+1-j
270  this%dof(k, 1, xh%lz, i) = edge_id
271  this%shared_dof(k, 1, xh%lz, i) = shared_dof
272  edge_id = edge_id + 1
273  end do
274 
275  call ep%edge_id(edge, 2)
276  shared_dof = msh%is_shared(edge)
277  global_id = msh%get_global(edge)
278  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
279  do j = 2, xh%lx - 1
280  k = j
281  if(int(edge%x(1), i8) .ne. this%dof(1,xh%ly,1,i)) k = xh%lx+1-j
282  this%dof(k, xh%ly, 1, i) = edge_id
283  this%shared_dof(k, xh%ly, 1, i) = shared_dof
284  edge_id = edge_id + 1
285  end do
286 
287  call ep%edge_id(edge, 4)
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  do j = 2, xh%lx - 1
292  k = j
293  if(int(edge%x(1), i8) .ne. this%dof(1,xh%ly,xh%lz,i)) k = xh%lx+1-j
294  this%dof(k, xh%ly, xh%lz, i) = edge_id
295  this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
296  edge_id = edge_id + 1
297  end do
298 
299 
300  !
301  ! Number edges in s-direction
302  !
303  call ep%edge_id(edge, 5)
304  shared_dof = msh%is_shared(edge)
305  global_id = msh%get_global(edge)
306  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
307  do j = 2, xh%ly - 1
308  k = j
309  if(int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) k = xh%ly+1-j
310  this%dof(1, k, 1, i) = edge_id
311  this%shared_dof(1, k, 1, i) = shared_dof
312  edge_id = edge_id + 1
313  end do
314 
315  call ep%edge_id(edge, 7)
316  shared_dof = msh%is_shared(edge)
317  global_id = msh%get_global(edge)
318  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
319  do j = 2, xh%ly - 1
320  k = j
321  if(int(edge%x(1), i8) .ne. this%dof(1,1,xh%lz,i)) k = xh%ly+1-j
322  this%dof(1, k, xh%lz, i) = edge_id
323  this%shared_dof(1, k, xh%lz, i) = shared_dof
324  edge_id = edge_id + 1
325  end do
326 
327  call ep%edge_id(edge, 6)
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  do j = 2, xh%ly - 1
332  k = j
333  if(int(edge%x(1), i8) .ne. this%dof(xh%lx,1,1,i)) k = xh%ly+1-j
334  this%dof(xh%lx, k, 1, i) = edge_id
335  this%shared_dof(xh%lx, k, 1, i) = shared_dof
336  edge_id = edge_id + 1
337  end do
338 
339  call ep%edge_id(edge, i8)
340  shared_dof = msh%is_shared(edge)
341  global_id = msh%get_global(edge)
342  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
343  do j = 2, xh%ly - 1
344  k = j
345  if(int(edge%x(1), i8) .ne. this%dof(xh%lx,1,xh%lz,i)) k = xh%lz+1-j
346  this%dof(xh%lx, k, xh%lz, i) = edge_id
347  this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
348  edge_id = edge_id + 1
349  end do
350 
351  !
352  ! Number edges in t-direction
353  !
354  call ep%edge_id(edge, 9)
355  shared_dof = msh%is_shared(edge)
356  global_id = msh%get_global(edge)
357  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
358  do j = 2, xh%lz - 1
359  k = j
360  if(int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) k = xh%lz+1-j
361  this%dof(1, 1, k, i) = edge_id
362  this%shared_dof(1, 1, k, i) = shared_dof
363  edge_id = edge_id + 1
364  end do
365 
366  call ep%edge_id(edge, 10)
367  shared_dof = msh%is_shared(edge)
368  global_id = msh%get_global(edge)
369  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
370  do j = 2, xh%lz - 1
371  k = j
372  if(int(edge%x(1), i8) .ne. this%dof(xh%lx,1,1,i)) k = xh%lz+1-j
373  this%dof(xh%lx, 1, k, i) = edge_id
374  this%shared_dof(xh%lx, 1, k, i) = shared_dof
375  edge_id = edge_id + 1
376  end do
377 
378  call ep%edge_id(edge, 11)
379  shared_dof = msh%is_shared(edge)
380  global_id = msh%get_global(edge)
381  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
382  do j = 2, xh%lz - 1
383  k = j
384  if(int(edge%x(1), i8) .ne. this%dof(1,xh%ly,1,i)) k = xh%lz+1-j
385  this%dof(1, xh%ly, k, i) = edge_id
386  this%shared_dof(1, xh%ly, k, i) = shared_dof
387  edge_id = edge_id + 1
388  end do
389 
390  call ep%edge_id(edge, 12)
391  shared_dof = msh%is_shared(edge)
392  global_id = msh%get_global(edge)
393  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(3)
394  do j = 2, xh%lz - 1
395  k = j
396  if(int(edge%x(1), i8) .ne. this%dof(xh%lx,xh%ly,1,i)) k = xh%lz+1-j
397  this%dof(xh%lx, xh%ly, k, i) = edge_id
398  this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
399  edge_id = edge_id + 1
400  end do
401  type is (quad_t)
402  !
403  ! Number edges in r-direction
404  !
405  call ep%facet_id(edge, 3)
406  shared_dof = msh%is_shared(edge)
407  global_id = msh%get_global(edge)
408  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
409  !Reverse order of tranversal if edge is reversed
410  do j = 2, xh%lx - 1
411  k = j
412  if(int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) k = xh%lx+1-j
413  this%dof(k, 1, 1, i) = edge_id
414  this%shared_dof(k, 1, 1, i) = shared_dof
415  edge_id = edge_id + 1
416  end do
417 
418  call ep%facet_id(edge, 4)
419  shared_dof = msh%is_shared(edge)
420  global_id = msh%get_global(edge)
421  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(1)
422  do j = 2, xh%lx - 1
423  k = j
424  if(int(edge%x(1), i8) .ne. this%dof(1,xh%ly,1,i)) k = xh%lx+1-j
425  this%dof(k, xh%ly, 1, i) = edge_id
426  this%shared_dof(k, xh%ly, 1, i) = shared_dof
427  edge_id = edge_id + 1
428  end do
429 
430  !
431  ! Number edges in s-direction
432  !
433  call ep%facet_id(edge, 1)
434  shared_dof = msh%is_shared(edge)
435  global_id = msh%get_global(edge)
436  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
437  do j = 2, xh%ly - 1
438  k = j
439  if(int(edge%x(1), i8) .ne. this%dof(1,1,1,i)) k = xh%ly+1-j
440  this%dof(1, k, 1, i) = edge_id
441  this%shared_dof(1, k, 1, i) = shared_dof
442  edge_id = edge_id + 1
443  end do
444 
445  call ep%facet_id(edge, 2)
446  shared_dof = msh%is_shared(edge)
447  global_id = msh%get_global(edge)
448  edge_id = edge_offset + int((global_id - 1), i8) * num_dofs_edges(2)
449  do j = 2, xh%ly - 1
450  k = j
451  if(int(edge%x(1), i8) .ne. this%dof(xh%lx,1,1,i)) k = xh%ly+1-j
452  this%dof(xh%lx, k, 1, i) = edge_id
453  this%shared_dof(xh%lx, k, 1, i) = shared_dof
454  edge_id = edge_id + 1
455  end do
456 
457  end select
458 
459  end do
460  end subroutine dofmap_number_edges
461 
463  subroutine dofmap_number_faces(this)
464  type(dofmap_t), target :: this
465  type(mesh_t), pointer :: msh
466  type(space_t), pointer :: Xh
467  integer :: i,j,k
468  integer :: global_id
469  type(tuple4_i4_t) :: face, face_order
470  integer(kind=i8) :: num_dofs_faces(3) ! #dofs for each dir (r, s, t)
471  integer(kind=i8) :: facet_offset, facet_id
472  logical :: shared_dof
473 
474  msh => this%msh
475  xh => this%Xh
476 
478  facet_offset = int(msh%glb_mpts, i8) + &
479  int(msh%glb_meds, i8) * int(xh%lx-2, i8) + int(1, i8)
480 
481  ! Number of dofs on an face excluding end-points
482  num_dofs_faces(1) = int((xh%ly - 2) * (xh%lz - 2), i8)
483  num_dofs_faces(2) = int((xh%lx - 2) * (xh%lz - 2), i8)
484  num_dofs_faces(3) = int((xh%lx - 2) * (xh%ly - 2), i8)
485 
486  do i = 1, msh%nelv
487 
488  !
489  ! Number facets in r-direction (s, t)-plane
490  !
491  call msh%elements(i)%e%facet_id(face, 1)
492  call msh%elements(i)%e%facet_order(face_order, 1)
493  shared_dof = msh%is_shared(face)
494  global_id = msh%get_global(face)
495  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(1)
496  do k = 2, xh%lz -1
497  do j = 2, xh%ly - 1
498  this%dof(1, j, k, i) = &
499  dofmap_facetidx(face_order,face,facet_id,j,k,xh%lz,xh%ly)
500  this%shared_dof(1, j, k, i) = shared_dof
501  end do
502  end do
503 
504  call msh%elements(i)%e%facet_id(face, 2)
505  call msh%elements(i)%e%facet_order(face_order, 2)
506  shared_dof = msh%is_shared(face)
507  global_id = msh%get_global(face)
508  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(1)
509  do k = 2, xh%lz -1
510  do j = 2, xh%ly - 1
511  this%dof(xh%lx, j, k, i) = &
512  dofmap_facetidx(face_order,face,facet_id,j,k,xh%lz,xh%ly)
513  this%shared_dof(xh%lx, j, k, i) = shared_dof
514  end do
515  end do
516 
517 
518  !
519  ! Number facets in s-direction (r, t)-plane
520  !
521  call msh%elements(i)%e%facet_id(face, 3)
522  call msh%elements(i)%e%facet_order(face_order, 3)
523  shared_dof = msh%is_shared(face)
524  global_id = msh%get_global(face)
525  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(2)
526  do k = 2, xh%lz - 1
527  do j = 2, xh%lx - 1
528  this%dof(j, 1, k, i) = &
529  dofmap_facetidx(face_order,face,facet_id,k,j,xh%lz,xh%lx)
530  this%shared_dof(j, 1, k, i) = shared_dof
531  end do
532  end do
533 
534  call msh%elements(i)%e%facet_id(face, 4)
535  call msh%elements(i)%e%facet_order(face_order, 4)
536  shared_dof = msh%is_shared(face)
537  global_id = msh%get_global(face)
538  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(2)
539  do k = 2, xh%lz - 1
540  do j = 2, xh%lx - 1
541  this%dof(j, xh%ly, k, i) = &
542  dofmap_facetidx(face_order,face,facet_id,k,j,xh%lz,xh%lx)
543  this%shared_dof(j, xh%ly, k, i) = shared_dof
544  end do
545  end do
546 
547 
548  !
549  ! Number facets in t-direction (r, s)-plane
550  !
551  call msh%elements(i)%e%facet_id(face, 5)
552  call msh%elements(i)%e%facet_order(face_order, 5)
553  shared_dof = msh%is_shared(face)
554  global_id = msh%get_global(face)
555  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(3)
556  do k = 2, xh%ly - 1
557  do j = 2, xh%lx - 1
558  this%dof(j, k, 1, i) = &
559  dofmap_facetidx(face_order,face,facet_id,k,j,xh%ly,xh%lx)
560  this%shared_dof(j, k, 1, i) = shared_dof
561  end do
562  end do
563 
564  call msh%elements(i)%e%facet_id(face, 6)
565  call msh%elements(i)%e%facet_order(face_order, 6)
566  shared_dof = msh%is_shared(face)
567  global_id = msh%get_global(face)
568  facet_id = facet_offset + int((global_id - 1), i8) * num_dofs_faces(3)
569  do k = 2, xh%ly - 1
570  do j = 2, xh%lx - 1
571  this%dof(j, k, xh%lz, i) = &
572  dofmap_facetidx(face_order,face,facet_id,k,j,xh%lz,xh%lx)
573  this%shared_dof(j, k, xh%lz, i) = shared_dof
574  end do
575  end do
576  end do
577 
578  end subroutine dofmap_number_faces
579 
581  function dofmap_facetidx(face_order, face, facet_id, k1, j1, lk1, lj1) result(facet_idx)
582  type(tuple4_i4_t) :: face_order, face
583  integer(kind=i8) :: facet_idx, facet_id
584  integer :: k1, j1, lk1, lj1
585  integer :: k,j,lk,lj
586 
587  k = k1 - 2
588  j = j1 - 2
589  lk = lk1 - 2
590  lj = lj1 - 2
591 
592  ! Given the indexes k,j for a GLL point on the inner part of the
593  ! face, we assign a unique number to it that depends on the
594  ! corner with the lowest id and its neighbour with the lowest
595  ! id. The id is assigned in this way to be consistent regardless
596  ! of how the faces are rotated or mirrored.
597  !
598  ! 4 -------- 3
599  ! | | k
600  ! |----->| ^
601  ! |----->| |
602  ! |----->| |
603  ! 1 -------- 2 0--->j
604 
605 
606  if(face_order%x(1) .eq. face%x(1)) then
607  if(face_order%x(2) .lt. face_order%x(4)) then
608  facet_idx = facet_id + j + k*lj
609  else
610  facet_idx = facet_id + j*lk + k
611  end if
612  else if(face_order%x(2) .eq. face%x(1)) then
613  if(face_order%x(3) .lt. face_order%x(1)) then
614  facet_idx = facet_id + lk*(lj-1-j) + k
615  else
616  facet_idx = facet_id + (lj-1-j) + k*lj
617  end if
618  else if(face_order%x(3) .eq. face%x(1)) then
619  if(face_order%x(4) .lt. face_order%x(2)) then
620  facet_idx = facet_id + (lj-1-j) + lj*(lk-1-k)
621  else
622  facet_idx = facet_id + lk*(lj-1-j) + (lk-1-k)
623  end if
624  else if(face_order%x(4) .eq. face%x(1)) then
625  if(face_order%x(1) .lt. face_order%x(3)) then
626  facet_idx = facet_id + lk*j + (lk-1-k)
627  else
628  facet_idx = facet_id + j + lj*(lk-1-k)
629  end if
630  end if
631 
632  end function dofmap_facetidx
633 
636  subroutine dofmap_generate_xyz(this)
637  type(dofmap_t), target :: this
638  integer :: i, j, el_idx
639  type(mesh_t), pointer :: msh
640  type(space_t), pointer :: Xh
641  real(kind=rp) :: rp_curve_data(5), curve_data_tot(5,12)
642  logical :: midpoint
643  integer :: n_edge, curve_type(12)
644 
645  msh => this%msh
646  xh => this%Xh
647 
648  if (msh%gdim .eq. 3) then
649  n_edge = 12
650  else
651  n_edge = 4
652  end if
653 
654  do i = 1, msh%nelv
655  call dofmap_xyzlin(xh, msh, msh%elements(i)%e, this%x(1,1,1,i), &
656  this%y(1,1,1,i), this%z(1,1,1,i))
657  end do
658  do i =1, msh%curve%size
659  midpoint = .false.
660  el_idx = msh%curve%curve_el(i)%el_idx
661  curve_type = msh%curve%curve_el(i)%curve_type
662  curve_data_tot = msh%curve%curve_el(i)%curve_data
663  do j = 1, n_edge
664  if (curve_type(j) .eq. 4) then
665  midpoint = .true.
666  end if
667  end do
668  if (midpoint .and. xh%lx .gt. 2) then
669  call dofmap_xyzquad(xh, msh, msh%elements(el_idx)%e, &
670  this%x(1,1,1,el_idx), this%y(1,1,1,el_idx),&
671  this%z(1,1,1,el_idx),curve_type, curve_data_tot)
672  end if
673  end do
674  do i =1, msh%curve%size
675  el_idx = msh%curve%curve_el(i)%el_idx
676  do j = 1, 8
677  if (msh%curve%curve_el(i)%curve_type(j) .eq. 3) then
678  rp_curve_data = msh%curve%curve_el(i)%curve_data(1:5,j)
679  call arc_surface(j, rp_curve_data, &
680  this%x(1,1,1,el_idx), &
681  this%y(1,1,1,el_idx), &
682  this%z(1,1,1, el_idx), &
683  xh, msh%elements(el_idx)%e, msh%gdim)
684  end if
685  end do
686  end do
687  if (associated(msh%apply_deform)) then
688  call msh%apply_deform(this%x, this%y, this%z, xh%lx, xh%ly, xh%lz)
689  end if
690  end subroutine dofmap_generate_xyz
691 
700  subroutine dofmap_xyzlin(Xh, msh, element, x, y, z)
701  type(mesh_t), pointer, intent(in) :: msh
702  type(space_t), intent(in) :: Xh
703  class(element_t), intent(in) :: element
704  real(kind=rp), intent(inout) :: x(xh%lx, xh%ly, xh%lz), &
705  y(xh%lx, xh%ly, xh%lz), &
706  z(xh%lx, xh%ly, xh%lz)
707  real(kind=rp) :: xyzb(2,2,2,3), zgml(xh%lx, 3)
708  real(kind=rp) :: jx(xh%lx*2)
709  real(kind=rp) :: jxt(xh%lx*2), jyt(xh%lx*2), jzt(xh%lx*2)
710  real(kind=rp) :: w(4*xh%lx**3), tmp(xh%lx, xh%lx, xh%lx)
711  real(kind=rp), dimension(2), parameter :: zlin = (/-1d0, 1d0/)
712 
713  integer :: j, k
714 
715  zgml = 0d0
716  xyzb = 0d0
717 
718  w = 0d0
719  call copy(zgml(1,1), xh%zg(1,1), xh%lx)
720  call copy(zgml(1,2), xh%zg(1,2), xh%ly)
721  if (msh%gdim .gt. 2) then
722  call copy(zgml(1,3), xh%zg(1,3), xh%lz)
723  end if
724 
725  k = 1
726  do j = 1, xh%lx
727  call fd_weights_full(zgml(j,1), zlin, 1, 0, jxt(k))
728  call fd_weights_full(zgml(j,2), zlin, 1, 0, jyt(k))
729  if (msh%gdim .gt. 2) then
730  call fd_weights_full(zgml(j,3), zlin, 1, 0, jzt(k))
731  end if
732  k = k + 2
733  end do
734  call trsp(jx, xh%lx, jxt, 2)
735 
736  if (msh%gdim .eq. 2) then
737  jzt = 1d0
738  end if
739 
740  do j = 1, msh%gdim
741  xyzb(1,1,1,j) = element%pts(1)%p%x(j)
742  xyzb(2,1,1,j) = element%pts(2)%p%x(j)
743  xyzb(1,2,1,j) = element%pts(3)%p%x(j)
744  xyzb(2,2,1,j) = element%pts(4)%p%x(j)
745 
746  if (msh%gdim .gt. 2) then
747  xyzb(1,1,2,j) = element%pts(5)%p%x(j)
748  xyzb(2,1,2,j) = element%pts(6)%p%x(j)
749  xyzb(1,2,2,j) = element%pts(7)%p%x(j)
750  xyzb(2,2,2,j) = element%pts(8)%p%x(j)
751  end if
752  end do
753  if (msh%gdim .eq. 3) then
754  call tensr3(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt, jzt, w)
755  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
756  call tensr3(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt, jzt, w)
757  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
758  call tensr3(tmp, xh%lz, xyzb(1,1,1,3), 2, jx, jyt, jzt, w)
759  call copy(z, tmp, xh%lx*xh%ly*xh%lz)
760  else
761  call tnsr2d_el(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt)
762  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
763  call tnsr2d_el(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt)
764  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
765  end if
766  end subroutine dofmap_xyzlin
767 
768  subroutine dofmap_xyzquad(Xh, msh, element, x, y, z, curve_type, curve_data)
769  type(mesh_t), pointer, intent(in) :: msh
770  type(space_t), intent(in) :: Xh
771  class(element_t), intent(in) :: element
772  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz), intent(inout) :: x, y, z
773  integer :: curve_type(12), eindx(12)
774  real(kind=rp) :: curve_data(5,12), x3(3,3,3), y3(3,3,3), z3(3,3,3)
775  type(space_t), target :: xh3
776  real(kind=rp), dimension(3), parameter :: zquad = (/-1d0, 0d0,1d0/)
777  real(kind=rp) :: zg(3)
778  real(kind=rp), dimension(Xh%lx,Xh%lx,Xh%lx) :: tmp
779  real(kind=rp) :: jx(xh%lx*3)
780  real(kind=rp) :: jxt(xh%lx*3), jyt(xh%lx*3), jzt(xh%lx*3)
781  real(kind=rp) :: w(4*xh%lxyz,2)
782  integer :: j, k, n_edges
783  eindx = (/2 , 6 , 8 , 4, &
784  20 , 24 , 26 , 22, &
785  10 , 12 , 18 , 16 /)
786 
787  w = 0d0
788  if (msh%gdim .eq. 3) then
789  n_edges = 12
790  call xh3%init(gll, 3, 3, 3)
791  else
792  n_edges = 4
793  call xh3%init(gll, 3, 3)
794  end if
795  call dofmap_xyzlin(xh3, msh, element, x3, y3, z3)
796 
797  do k = 1, n_edges
798  if (curve_type(k) .eq. 4) then
799  x3(eindx(k),1,1) = curve_data(1,k)
800  y3(eindx(k),1,1) = curve_data(2,k)
801  z3(eindx(k),1,1) = curve_data(3,k)
802  end if
803  end do
804  zg(1) = -1
805  zg(2) = 0
806  zg(3) = 1
807  if (msh%gdim .eq. 3) then
808  call gh_face_extend_3d(x3,zg,3,2,w(1,1),w(1,2)) ! 2 --> edge extend
809  call gh_face_extend_3d(y3,zg,3,2,w(1,1),w(1,2))
810  call gh_face_extend_3d(z3,zg,3,2,w(1,1),w(1,2))
811  else
812  call neko_warning(' m deformation not supported for 2d yet')
813  call gh_face_extend_2d(x3,zg,3,2,w(1,1),w(1,2)) ! 2 --> edge extend
814  call gh_face_extend_2d(y3,zg,3,2,w(1,1),w(1,2))
815  end if
816  k =1
817  do j = 1, xh%lx
818  call fd_weights_full(xh%zg(j,1),zquad,2,0,jxt(k))
819  call fd_weights_full(xh%zg(j,2),zquad,2,0,jyt(k))
820  if (msh%gdim .gt. 2) then
821  call fd_weights_full(xh%zg(j,3),zquad,2,0,jzt(k))
822  end if
823  k = k + 3
824  end do
825  call trsp(jx, xh%lx, jxt, 3)
826  if (msh%gdim .eq. 3) then
827  call tensr3(tmp, xh%lx, x3, 3, jx, jyt, jzt, w)
828  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
829  call tensr3(tmp, xh%ly, y3, 3, jx, jyt, jzt, w)
830  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
831  call tensr3(tmp, xh%lz, z3, 3, jx, jyt, jzt, w)
832  call copy(z, tmp, xh%lx*xh%ly*xh%lz)
833  else
834  call tnsr2d_el(tmp, xh%lx, x3, 3, jx, jyt)
835  call copy(x, tmp, xh%lx*xh%ly*xh%lz)
836  call tnsr2d_el(tmp, xh%ly, y3, 3, jx, jyt)
837  call copy(y, tmp, xh%lx*xh%ly*xh%lz)
838  end if
839 
840  call xh3%free()
841  end subroutine dofmap_xyzquad
842 
848  subroutine gh_face_extend_3d(x, zg, n, gh_type, e, v)
849  integer, intent(in) :: n
850  real(kind=rp), intent(inout) :: x(n, n, n)
851  real(kind=rp), intent(in) :: zg(n)
852  real(kind=rp), intent(inout) :: e(n, n, n)
853  real(kind=rp), intent(inout) :: v(n, n, n)
854  integer :: gh_type, ntot, kk, jj, ii, k, j, i
855  real(kind=rp) :: si, sj, sk, hi, hj, hk
856 
857 !
858 ! Build vertex interpolant
859 !
860  ntot = n**3
861  call rzero(v, ntot)
862  do kk = 1, n, n-1
863  do jj = 1, n, n-1
864  do ii = 1, n, n-1
865  do k = 1, n
866  do j = 1, n
867  do i = 1, n
868  si = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
869  sj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
870  sk = 0.5*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
871  v(i,j,k) = v(i,j,k) + si*sj*sk*x(ii,jj,kk)
872  end do
873  end do
874  end do
875  end do
876  end do
877  end do
878 
879  if (gh_type .eq. 1) then
880  call copy(x, v, ntot)
881  return
882  end if
883 !
884 !
885 ! Extend 12 edges
886  call rzero(e, ntot)
887 !
888 ! x-edges
889 !
890  do kk = 1, n, n-1
891  do jj = 1, n, n-1
892  do k = 1, n
893  do j = 1, n
894  do i = 1, n
895  hj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
896  hk = 0.5*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
897  e(i,j,k) = e(i,j,k) + hj*hk*(x(i,jj,kk)-v(i,jj,kk))
898  end do
899  end do
900  end do
901  end do
902  end do
903 !
904 ! y-edges
905 !
906  do kk = 1, n, n-1
907  do ii = 1, n, n-1
908  do k = 1, n
909  do j = 1, n
910  do i = 1, n
911  hi = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
912  hk = 0.5*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
913  e(i,j,k) = e(i,j,k) + hi*hk*(x(ii,j,kk)-v(ii,j,kk))
914  end do
915  end do
916  end do
917  end do
918  end do
919 !
920 ! z-edges
921 !
922  do jj = 1, n, n-1
923  do ii = 1, n, n-1
924  do k = 1, n
925  do j = 1, n
926  do i = 1, n
927  hi = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
928  hj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
929  e(i,j,k) = e(i,j,k) + hi*hj*(x(ii,jj,k)-v(ii,jj,k))
930  end do
931  end do
932  end do
933  end do
934  end do
935 
936  call add2(e, v, ntot)
937 
938  if (gh_type .eq. 2) then
939  call copy(x, e, ntot)
940  return
941  end if
942 !
943 ! Extend faces
944 !
945  call rzero(v, ntot)
946 !
947 ! x-edges
948 !
949  do ii = 1, n, n-1
950  do k = 1, n
951  do j = 1, n
952  do i = 1, n
953  hi = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
954  v(i,j,k) = v(i,j,k) + hi*(x(ii,j,k)-e(ii,j,k))
955  end do
956  end do
957  end do
958  end do
959 !
960 ! y-edges
961 !
962  do jj = 1, n, n-1
963  do k = 1, n
964  do j = 1, n
965  do i = 1, n
966  hj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
967  v(i,j,k) = v(i,j,k) + hj*(x(i,jj,k)-e(i,jj,k))
968  end do
969  end do
970  end do
971  end do
972 !
973 ! z-edges
974 !
975  do kk= 1 , n, n-1
976  do k = 1, n
977  do j = 1, n
978  do i = 1, n
979  hk = 0.5*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
980  v(i,j,k) = v(i,j,k) + hk*(x(i,j,kk)-e(i,j,kk))
981  end do
982  end do
983  end do
984  end do
985 
986  call add2(v, e, ntot)
987  call copy(x, v ,ntot)
988 
989  end subroutine gh_face_extend_3d
990 
994  subroutine gh_face_extend_2d(x, zg, n, gh_type, e, v)
995  integer, intent(in) :: n
996  real(kind=rp), intent(inout) :: x(n, n)
997  real(kind=rp), intent(in) :: zg(n)
998  real(kind=rp), intent(inout) :: e(n, n)
999  real(kind=rp), intent(inout) :: v(n, n)
1000  integer, intent(in) :: gh_type
1001  integer :: i,j , jj, ii, ntot
1002  real(kind=rp) :: si, sj, hi, hj
1003 
1004  !Build vertex interpolant
1005 
1006  ntot=n*n
1007  call rzero(v,ntot)
1008  do jj = 1, n, n-1
1009  do ii = 1, n, n-1
1010  do j = 1, n
1011  do i = 1, n
1012  si = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1013  sj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1014  v(i,j) = v(i,j) + si*sj*x(ii,jj)
1015  end do
1016  end do
1017  end do
1018  end do
1019  if (gh_type .eq. 1) then
1020  call copy(x,v,ntot)
1021  return
1022  end if
1023 
1024  !Extend 4 edges
1025  call rzero(e,ntot)
1026 
1027  !x-edges
1028 
1029  do jj = 1, n, n-1
1030  do j = 1, n
1031  do i = 1, n
1032  hj = 0.5*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1033  e(i,j) = e(i,j) + hj*(x(i,jj)-v(i,jj))
1034  end do
1035  end do
1036  end do
1037 
1038  !y-edges
1039 
1040  do ii = 1, n, n-1
1041  do j = 1, n
1042  do i = 1, n
1043  hi = 0.5*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1044  e(i,j) = e(i,j) + hi*(x(ii,j)-v(ii,j))
1045  end do
1046  end do
1047  end do
1048 
1049  call add3(x, e, v, ntot)
1050 
1051  end subroutine gh_face_extend_2d
1052 
1053 
1054 
1055  subroutine arc_surface(isid, curve_data, x, y, z, Xh, element, gdim)
1056  integer, intent(in) :: isid, gdim
1057  type(space_t), intent(in) :: Xh
1058  class(element_t) :: element
1059  real(kind=rp), dimension(5), intent(in) :: curve_data
1060  real(kind=rp), dimension(Xh%lx, Xh%ly, Xh%lz), intent(inout) :: x, y, z
1061  real(kind=rp) :: pt1x, pt1y, pt2x, pt2y, pt12x, pt12y
1062  real(kind=rp) :: radius, dtheta, r, xys
1063  real(kind=rp) :: theta0, xcenn, ycenn, h(xh%lx, 3, 2)
1064  real(kind=rp) :: xcrved(xh%lx), ycrved(xh%lx), xs, ys
1065  integer :: isid1, ixt, iyt, izt, ix, itmp
1066  integer(i4), dimension(6), parameter :: fcyc_to_sym = (/3, 2, 4, 1, 5, 6/) ! cyclic to symmetric face mapping
1067  integer(i4), dimension(12), parameter :: ecyc_to_sym = (/1, 6, 2, 5, 3, 8, &
1068  & 4, 7, 9, 10, 12, 11/) ! cyclic to symmetric edge mapping
1069  integer, parameter, dimension(2, 12) :: edge_nodes = reshape((/1, 2, 3, 4, 5, 6, &
1070  & 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
1071  ! copy from hex as this has private attribute there
1072 
1073  ! this subroutine is a mess of symmetric and cyclic edge/face numberring and
1074  ! cannot be cleaned without changing an input format (isid seems to be
1075  ! a cyclic edge number)
1076  ! following according to cyclic edge numbering and orientation
1077  itmp = ecyc_to_sym(isid)
1078  select case(isid)
1079  case(1:2,5:6)
1080  pt1x = element%pts(edge_nodes(1,itmp))%p%x(1)
1081  pt1y = element%pts(edge_nodes(1,itmp))%p%x(2)
1082  pt2x = element%pts(edge_nodes(2,itmp))%p%x(1)
1083  pt2y = element%pts(edge_nodes(2,itmp))%p%x(2)
1084  case(3:4,7:8)
1085  pt1x = element%pts(edge_nodes(2,itmp))%p%x(1)
1086  pt1y = element%pts(edge_nodes(2,itmp))%p%x(2)
1087  pt2x = element%pts(edge_nodes(1,itmp))%p%x(1)
1088  pt2y = element%pts(edge_nodes(1,itmp))%p%x(2)
1089  end select
1090  ! find slope of perpendicular
1091  radius = curve_data(1)
1092  xs = pt2y-pt1y
1093  ys = pt1x-pt2x
1094  ! make length radius
1095  xys = sqrt(xs**2 + ys**2)
1096  ! sanity check
1097  if (abs(2.0 * radius) <= xys * 1.00001) &
1098  & call neko_error('Radius to small for arced element surface')
1099  ! find center
1100  dtheta = abs(asin(0.5*xys/radius))
1101  pt12x = (pt1x + pt2x)/2.0
1102  pt12y = (pt1y + pt2y)/2.0
1103  xcenn = pt12x - xs/xys * radius*cos(dtheta)
1104  ycenn = pt12y - ys/xys * radius*cos(dtheta)
1105  theta0 = atan2((pt12y-ycenn), (pt12x-xcenn))
1106 ! compute perturbation of geometry
1107  isid1 = mod(isid+4-1, 4)+1
1108  call compute_h(h, xh%zg, gdim, xh%lx)
1109  if (radius < 0.0) dtheta = -dtheta
1110  do ix=1,xh%lx
1111  ixt=ix
1112  if (isid1.gt.2) ixt=xh%lx+1-ix
1113  r=xh%zg(ix,1)
1114  xcrved(ixt) = xcenn + abs(radius) * cos(theta0 + r*dtheta) &
1115  - ( h(ix,1,1)*pt1x + h(ix,1,2)*pt2x )
1116  ycrved(ixt) = ycenn + abs(radius) * sin(theta0 + r*dtheta) &
1117  - ( h(ix,1,1)*pt1y + h(ix,1,2)*pt2y )
1118  end do
1119 ! points all set, add perturbation to current mesh.
1120 ! LEGACY WARNING
1121 ! I dont want to dive in this again, Martin Karp 2/3 - 2021
1122  isid1 = fcyc_to_sym(isid1)
1123  izt = (isid-1)/4+1
1124  iyt = isid1-2
1125  ixt = isid1
1126  if (isid1 .le. 2) then
1127  call addtnsr(x, h(1,1,ixt), xcrved, h(1,3,izt) &
1128  ,xh%lx, xh%ly, xh%lz)
1129  call addtnsr(y, h(1,1,ixt), ycrved, h(1,3,izt) &
1130  ,xh%lx, xh%ly, xh%lz)
1131  else
1132  call addtnsr(x, xcrved, h(1,2,iyt), h(1,3,izt) &
1133  ,xh%lx, xh%ly, xh%lz)
1134  call addtnsr(y, ycrved, h(1,2,iyt), h(1,3,izt) &
1135  ,xh%lx, xh%ly, xh%lz)
1136  end if
1137  end subroutine arc_surface
1138 
1139  subroutine compute_h(h, zgml, gdim, lx)
1140  integer, intent(in) :: lx, gdim
1141  real(kind=rp), intent(inout) :: h(lx, 3, 2)
1142  real(kind=rp), intent(in) :: zgml(lx, 3)
1143  integer :: ix, iy, iz
1144 
1145  do ix = 1, lx
1146  h(ix,1,1) = (1.0_rp - zgml(ix, 1)) * 0.5_rp
1147  h(ix,1,2) = (1.0_rp + zgml(ix, 1)) * 0.5_rp
1148  end do
1149 
1150  do iy = 1, lx
1151  h(iy,2,1) = (1.0_rp - zgml(iy, 2)) * 0.5_rp
1152  h(iy,2,2) = (1.0_rp + zgml(iy, 2)) * 0.5_rp
1153  end do
1154 
1155  if (gdim .eq. 3) then
1156  do iz = 1, lx
1157  h(iz,3,1) = (1.0_rp - zgml(iz, 3)) * 0.5_rp
1158  h(iz,3,2) = (1.0_rp + zgml(iz, 3)) * 0.5_rp
1159  end do
1160  else
1161  call rone(h(1,3,1), lx)
1162  call rone(h(1,3,2), lx)
1163  end if
1164 
1165  end subroutine compute_h
1166 
1167 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:172
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:995
subroutine dofmap_xyzquad(Xh, msh, element, x, y, z, curve_type, curve_data)
Definition: dofmap.f90:769
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:637
subroutine compute_h(h, zgml, gdim, lx)
Definition: dofmap.f90:1140
subroutine arc_surface(isid, curve_data, x, y, z, Xh, element, gdim)
Definition: dofmap.f90:1056
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
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:582
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:701
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:464
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:849
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:200
subroutine, public add3(a, b, c, n)
Vector addition .
Definition: math.f90:516
subroutine, public add2(a, b, n)
Vector addition .
Definition: math.f90:503
subroutine, public copy(a, b, n)
Copy a vector .
Definition: math.f90:211
subroutine, public rzero(a, n)
Zero a real vector.
Definition: math.f90:167
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 neko_warning(warning_msg)
Definition: utils.f90:191
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