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