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