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