49 use,
intrinsic :: iso_c_binding, only : c_ptr, c_null_ptr
54 integer(kind=i8),
allocatable :: dof(:,:,:,:)
55 logical,
allocatable :: shared_dof(:,:,:,:)
56 real(kind=
rp),
allocatable :: x(:,:,:,:)
57 real(kind=
rp),
allocatable :: y(:,:,:,:)
58 real(kind=
rp),
allocatable :: z(:,:,:,:)
59 integer,
private :: ntot
67 type(c_ptr) :: x_d = c_null_ptr
68 type(c_ptr) :: y_d = c_null_ptr
69 type(c_ptr) :: z_d = c_null_ptr
87 type(
mesh_t),
target,
intent(inout) :: msh
88 type(
space_t),
target,
intent(inout) :: Xh
90 if ((msh%gdim .eq. 3 .and. xh%lz .eq. 1) .or. &
91 (msh%gdim .eq. 2 .and. xh%lz .gt. 1))
then
92 call neko_error(
"Invalid dimension of function space for the given mesh")
100 this%ntot = xh%lx* xh%ly * xh%lz * msh%nelv
106 allocate(this%dof(xh%lx, xh%ly, xh%lz, msh%nelv))
107 allocate(this%shared_dof(xh%lx, xh%ly, xh%lz, msh%nelv))
110 this%shared_dof = .false.
113 if (msh%gdim .eq. 3)
then
126 allocate(this%x(xh%lx, xh%ly, xh%lz, msh%nelv))
127 allocate(this%y(xh%lx, xh%ly, xh%lz, msh%nelv))
128 allocate(this%z(xh%lx, xh%ly, xh%lz, msh%nelv))
154 class(
dofmap_t),
intent(inout) :: this
156 if (
allocated(this%dof))
then
160 if (
allocated(this%shared_dof))
then
161 deallocate(this%shared_dof)
164 if (
allocated(this%x))
then
168 if (
allocated(this%y))
then
172 if (
allocated(this%z))
then
182 if (c_associated(this%x_d))
then
186 if (c_associated(this%y_d))
then
190 if (c_associated(this%z_d))
then
206 integer :: il, jl, ix, iy, iz
207 type(
mesh_t),
pointer :: msh
214 ix = mod(jl - 1, 2) * (xh%lx - 1) + 1
215 iy = (mod(jl - 1, 4)/2) * (xh%ly - 1) + 1
216 iz = ((jl - 1)/4) * (xh%lz - 1) + 1
217 this%dof(ix, iy, iz, il) = int(msh%elements(il)%e%pts(jl)%p%id(),
i8)
218 this%shared_dof(ix, iy, iz, il) = &
219 msh%is_shared(msh%elements(il)%e%pts(jl)%p)
227 type(
mesh_t),
pointer :: msh
232 integer(kind=i8) :: num_dofs_edges(3)
233 integer(kind=i8) :: edge_id, edge_offset
234 logical :: shared_dof
240 num_dofs_edges(1) = int(xh%lx - 2,
i8)
241 num_dofs_edges(2) = int(xh%ly - 2,
i8)
242 num_dofs_edges(3) = int(xh%lz - 2,
i8)
243 edge_offset = int(msh%glb_mpts,
i8) + int(1, 4)
247 select type (ep => msh%elements(i)%e)
252 call ep%edge_id(edge, 1)
253 shared_dof = msh%is_shared(edge)
254 global_id = msh%get_global(edge)
255 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
257 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
258 do concurrent(j = 2:xh%lx - 1)
260 this%dof(k, 1, 1, i) = edge_id + (j-2)
261 this%shared_dof(k, 1, 1, i) = shared_dof
264 do concurrent(j = 2:xh%lx - 1)
266 this%dof(k, 1, 1, i) = edge_id + (j-2)
267 this%shared_dof(k, 1, 1, i) = shared_dof
271 call ep%edge_id(edge, 3)
272 shared_dof = msh%is_shared(edge)
273 global_id = msh%get_global(edge)
274 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
275 if (int(edge%x(1),
i8) .ne. this%dof(1, 1, xh%lz, i))
then
276 do concurrent(j = 2:xh%lx - 1)
278 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
279 this%shared_dof(k, 1, xh%lz, i) = shared_dof
282 do concurrent(j = 2:xh%lx - 1)
284 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
285 this%shared_dof(k, 1, xh%lz, i) = shared_dof
289 call ep%edge_id(edge, 2)
290 shared_dof = msh%is_shared(edge)
291 global_id = msh%get_global(edge)
292 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
293 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
294 do concurrent(j = 2:xh%lx - 1)
296 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
297 this%shared_dof(k, xh%ly, 1, i) = shared_dof
300 do concurrent(j = 2:xh%lx - 1)
302 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
303 this%shared_dof(k, xh%ly, 1, i) = shared_dof
307 call ep%edge_id(edge, 4)
308 shared_dof = msh%is_shared(edge)
309 global_id = msh%get_global(edge)
310 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
311 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, xh%lz, i))
then
312 do concurrent(j = 2:xh%lx - 1)
314 this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
315 this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
318 do concurrent(j = 2:xh%lx - 1)
320 this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
321 this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
329 call ep%edge_id(edge, 5)
330 shared_dof = msh%is_shared(edge)
331 global_id = msh%get_global(edge)
332 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
333 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
334 do concurrent(j = 2:xh%ly - 1)
336 this%dof(1, k, 1, i) = edge_id + (j-2)
337 this%shared_dof(1, k, 1, i) = shared_dof
340 do concurrent(j = 2:xh%ly - 1)
342 this%dof(1, k, 1, i) = edge_id + (j-2)
343 this%shared_dof(1, k, 1, i) = shared_dof
347 call ep%edge_id(edge, 7)
348 shared_dof = msh%is_shared(edge)
349 global_id = msh%get_global(edge)
350 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
351 if (int(edge%x(1),
i8) .ne. this%dof(1, 1, xh%lz, i))
then
352 do concurrent(j = 2:xh%ly - 1)
354 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
355 this%shared_dof(1, k, xh%lz, i) = shared_dof
358 do concurrent(j = 2:xh%ly - 1)
360 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
361 this%shared_dof(1, k, xh%lz, i) = shared_dof
365 call ep%edge_id(edge, 6)
366 shared_dof = msh%is_shared(edge)
367 global_id = msh%get_global(edge)
368 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
369 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, 1, 1, i))
then
370 do concurrent(j = 2:xh%ly - 1)
372 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
373 this%shared_dof(xh%lx, k, 1, i) = shared_dof
376 do concurrent(j = 2:xh%ly - 1)
378 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
379 this%shared_dof(xh%lx, k, 1, i) = shared_dof
383 call ep%edge_id(edge,
i8)
384 shared_dof = msh%is_shared(edge)
385 global_id = msh%get_global(edge)
386 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
387 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, 1, xh%lz, i))
then
388 do concurrent(j = 2:xh%ly - 1)
390 this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
391 this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
394 do concurrent(j = 2:xh%ly - 1)
396 this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
397 this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
404 call ep%edge_id(edge, 9)
405 shared_dof = msh%is_shared(edge)
406 global_id = msh%get_global(edge)
407 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
408 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
409 do concurrent(j = 2:xh%lz - 1)
411 this%dof(1, 1, k, i) = edge_id + (j-2)
412 this%shared_dof(1, 1, k, i) = shared_dof
415 do concurrent(j = 2:xh%lz - 1)
417 this%dof(1, 1, k, i) = edge_id + (j-2)
418 this%shared_dof(1, 1, k, i) = shared_dof
422 call ep%edge_id(edge, 10)
423 shared_dof = msh%is_shared(edge)
424 global_id = msh%get_global(edge)
425 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
426 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx,1,1,i))
then
427 do concurrent(j = 2:xh%lz - 1)
429 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
430 this%shared_dof(xh%lx, 1, k, i) = shared_dof
433 do concurrent(j = 2:xh%lz - 1)
435 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
436 this%shared_dof(xh%lx, 1, k, i) = shared_dof
440 call ep%edge_id(edge, 11)
441 shared_dof = msh%is_shared(edge)
442 global_id = msh%get_global(edge)
443 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
444 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
445 do concurrent(j = 2:xh%lz - 1)
447 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
448 this%shared_dof(1, xh%ly, k, i) = shared_dof
451 do concurrent(j = 2:xh%lz - 1)
453 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
454 this%shared_dof(1, xh%ly, k, i) = shared_dof
458 call ep%edge_id(edge, 12)
459 shared_dof = msh%is_shared(edge)
460 global_id = msh%get_global(edge)
461 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
462 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, xh%ly, 1, i))
then
463 do concurrent(j = 2:xh%lz - 1)
465 this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
466 this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
469 do concurrent(j = 2:xh%lz - 1)
471 this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
472 this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
479 call ep%facet_id(edge, 3)
480 shared_dof = msh%is_shared(edge)
481 global_id = msh%get_global(edge)
482 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
484 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
485 do concurrent(j = 2:xh%lx - 1)
487 this%dof(k, 1, 1, i) = edge_id + (j-2)
488 this%shared_dof(k, 1, 1, i) = shared_dof
491 do concurrent(j = 2:xh%lx - 1)
493 this%dof(k, 1, 1, i) = edge_id + (j-2)
494 this%shared_dof(k, 1, 1, i) = shared_dof
498 call ep%facet_id(edge, 4)
499 shared_dof = msh%is_shared(edge)
500 global_id = msh%get_global(edge)
501 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
502 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
503 do concurrent(j = 2:xh%lx - 1)
505 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
506 this%shared_dof(k, xh%ly, 1, i) = shared_dof
509 do concurrent(j = 2:xh%lx - 1)
511 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
512 this%shared_dof(k, xh%ly, 1, i) = shared_dof
519 call ep%facet_id(edge, 1)
520 shared_dof = msh%is_shared(edge)
521 global_id = msh%get_global(edge)
522 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
523 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
524 do concurrent(j = 2:xh%ly - 1)
526 this%dof(1, k, 1, i) = edge_id + (j-2)
527 this%shared_dof(1, k, 1, i) = shared_dof
530 do concurrent(j = 2:xh%ly - 1)
532 this%dof(1, k, 1, i) = edge_id + (j-2)
533 this%shared_dof(1, k, 1, i) = shared_dof
537 call ep%facet_id(edge, 2)
538 shared_dof = msh%is_shared(edge)
539 global_id = msh%get_global(edge)
540 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
541 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx,1,1,i))
then
542 do concurrent(j = 2:xh%ly - 1)
544 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
545 this%shared_dof(xh%lx, k, 1, i) = shared_dof
548 do concurrent(j = 2:xh%ly - 1)
550 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
551 this%shared_dof(xh%lx, k, 1, i) = shared_dof
562 type(
mesh_t),
pointer :: msh
567 integer(kind=i8) :: num_dofs_faces(3)
568 integer(kind=i8) :: facet_offset, facet_id
569 logical :: shared_dof
575 facet_offset = int(msh%glb_mpts,
i8) + &
576 int(msh%glb_meds,
i8) * int(xh%lx-2,
i8) + int(1,
i8)
579 num_dofs_faces(1) = int((xh%ly - 2) * (xh%lz - 2),
i8)
580 num_dofs_faces(2) = int((xh%lx - 2) * (xh%lz - 2),
i8)
581 num_dofs_faces(3) = int((xh%lx - 2) * (xh%ly - 2),
i8)
588 call msh%elements(i)%e%facet_id(face, 1)
589 call msh%elements(i)%e%facet_order(face_order, 1)
590 shared_dof = msh%is_shared(face)
591 global_id = msh%get_global(face)
592 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(1)
593 do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
594 this%dof(1, j, k, i) = &
596 this%shared_dof(1, j, k, i) = shared_dof
599 call msh%elements(i)%e%facet_id(face, 2)
600 call msh%elements(i)%e%facet_order(face_order, 2)
601 shared_dof = msh%is_shared(face)
602 global_id = msh%get_global(face)
603 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(1)
604 do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
605 this%dof(xh%lx, j, k, i) = &
607 this%shared_dof(xh%lx, j, k, i) = shared_dof
614 call msh%elements(i)%e%facet_id(face, 3)
615 call msh%elements(i)%e%facet_order(face_order, 3)
616 shared_dof = msh%is_shared(face)
617 global_id = msh%get_global(face)
618 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(2)
619 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
620 this%dof(j, 1, k, i) = &
622 this%shared_dof(j, 1, k, i) = shared_dof
625 call msh%elements(i)%e%facet_id(face, 4)
626 call msh%elements(i)%e%facet_order(face_order, 4)
627 shared_dof = msh%is_shared(face)
628 global_id = msh%get_global(face)
629 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(2)
630 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
631 this%dof(j, xh%ly, k, i) = &
633 this%shared_dof(j, xh%ly, k, i) = shared_dof
640 call msh%elements(i)%e%facet_id(face, 5)
641 call msh%elements(i)%e%facet_order(face_order, 5)
642 shared_dof = msh%is_shared(face)
643 global_id = msh%get_global(face)
644 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(3)
645 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
646 this%dof(j, k, 1, i) = &
648 this%shared_dof(j, k, 1, i) = shared_dof
651 call msh%elements(i)%e%facet_id(face, 6)
652 call msh%elements(i)%e%facet_order(face_order, 6)
653 shared_dof = msh%is_shared(face)
654 global_id = msh%get_global(face)
655 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(3)
656 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
657 this%dof(j, k, xh%lz, i) = &
659 this%shared_dof(j, k, xh%lz, i) = shared_dof
667 lj1)
result(facet_idx)
669 integer(kind=i8),
intent(in) :: facet_id
670 integer(kind=i8) :: facet_idx
671 integer,
intent(in) :: k1, j1, lk1, lj1
672 integer :: k, j, lk, lj
693 if (face_order%x(1) .eq. face%x(1))
then
694 if (face_order%x(2) .lt. face_order%x(4))
then
695 facet_idx = facet_id + j + k*lj
697 facet_idx = facet_id + j*lk + k
699 else if (face_order%x(2) .eq. face%x(1))
then
700 if (face_order%x(3) .lt. face_order%x(1))
then
701 facet_idx = facet_id + lk*(lj-1-j) + k
703 facet_idx = facet_id + (lj-1-j) + k*lj
705 else if (face_order%x(3) .eq. face%x(1))
then
706 if (face_order%x(4) .lt. face_order%x(2))
then
707 facet_idx = facet_id + (lj-1-j) + lj*(lk-1-k)
709 facet_idx = facet_id + lk*(lj-1-j) + (lk-1-k)
711 else if (face_order%x(4) .eq. face%x(1))
then
712 if (face_order%x(1) .lt. face_order%x(3))
then
713 facet_idx = facet_id + lk*j + (lk-1-k)
715 facet_idx = facet_id + j + lj*(lk-1-k)
725 integer :: i, j, el_idx
726 type(
mesh_t),
pointer :: msh
728 real(kind=
rp) :: rp_curve_data(5), curve_data_tot(5,12)
730 integer :: n_edge, curve_type(12)
735 if (msh%gdim .eq. 3)
then
742 call dofmap_xyzlin(xh, msh, msh%elements(i)%e, this%x(1,1,1,i), &
743 this%y(1,1,1,i), this%z(1,1,1,i))
745 do i = 1, msh%curve%size
747 el_idx = msh%curve%curve_el(i)%el_idx
748 curve_type = msh%curve%curve_el(i)%curve_type
749 curve_data_tot = msh%curve%curve_el(i)%curve_data
751 if (curve_type(j) .eq. 4)
then
755 if (midpoint .and. xh%lx .gt. 2)
then
757 this%x(1, 1, 1, el_idx), this%y(1, 1, 1, el_idx), &
758 this%z(1 ,1, 1, el_idx), curve_type, curve_data_tot)
761 do i = 1, msh%curve%size
762 el_idx = msh%curve%curve_el(i)%el_idx
764 if (msh%curve%curve_el(i)%curve_type(j) .eq. 3)
then
765 rp_curve_data = msh%curve%curve_el(i)%curve_data(1:5,j)
767 this%x(1, 1, 1, el_idx), &
768 this%y(1, 1, 1, el_idx), &
769 this%z(1, 1, 1, el_idx), &
770 xh, msh%elements(el_idx)%e, msh%gdim)
774 if (
associated(msh%apply_deform))
then
775 call msh%apply_deform(this%x, this%y, this%z, xh%lx, xh%ly, xh%lz)
788 type(
mesh_t),
pointer,
intent(in) :: msh
789 type(
space_t),
intent(in) :: Xh
791 real(kind=
rp),
intent(inout) :: x(xh%lx, xh%ly, xh%lz), &
792 y(xh%lx, xh%ly, xh%lz), &
793 z(xh%lx, xh%ly, xh%lz)
794 real(kind=
rp) :: xyzb(2,2,2,3), zgml(xh%lx, 3)
795 real(kind=
rp) :: jx(xh%lx*2)
796 real(kind=
rp) :: jxt(xh%lx*2), jyt(xh%lx*2), jzt(xh%lx*2)
797 real(kind=
rp) :: w(4*xh%lx**3), tmp(xh%lx, xh%lx, xh%lx)
798 real(kind=
rp),
dimension(2),
parameter :: zlin = [-1d0, 1d0]
806 call copy(zgml(1,1), xh%zg(1,1), xh%lx)
807 call copy(zgml(1,2), xh%zg(1,2), xh%ly)
808 if (msh%gdim .gt. 2)
then
809 call copy(zgml(1,3), xh%zg(1,3), xh%lz)
816 if (msh%gdim .gt. 2)
then
821 call trsp(jx, xh%lx, jxt, 2)
823 if (msh%gdim .eq. 2)
then
827 if (msh%gdim .gt. 2)
then
828 do concurrent(j = 1:msh%gdim)
829 xyzb(1,1,1,j) =
element%pts(1)%p%x(j)
830 xyzb(2,1,1,j) =
element%pts(2)%p%x(j)
831 xyzb(1,2,1,j) =
element%pts(3)%p%x(j)
832 xyzb(2,2,1,j) =
element%pts(4)%p%x(j)
834 xyzb(1,1,2,j) =
element%pts(5)%p%x(j)
835 xyzb(2,1,2,j) =
element%pts(6)%p%x(j)
836 xyzb(1,2,2,j) =
element%pts(7)%p%x(j)
837 xyzb(2,2,2,j) =
element%pts(8)%p%x(j)
840 do concurrent(j = 1:msh%gdim)
841 xyzb(1,1,1,j) =
element%pts(1)%p%x(j)
842 xyzb(2,1,1,j) =
element%pts(2)%p%x(j)
843 xyzb(1,2,1,j) =
element%pts(3)%p%x(j)
844 xyzb(2,2,1,j) =
element%pts(4)%p%x(j)
847 if (msh%gdim .eq. 3)
then
848 call tensr3(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt, jzt, w)
849 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
850 call tensr3(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt, jzt, w)
851 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
852 call tensr3(tmp, xh%lz, xyzb(1,1,1,3), 2, jx, jyt, jzt, w)
853 call copy(z, tmp, xh%lx*xh%ly*xh%lz)
855 call tnsr2d_el(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt)
856 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
857 call tnsr2d_el(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt)
858 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
863 type(
mesh_t),
pointer,
intent(in) :: msh
864 type(
space_t),
intent(in) :: Xh
866 real(kind=
rp),
dimension(Xh%lx, Xh%ly, Xh%lz),
intent(inout) :: x, y, z
867 integer :: curve_type(12), eindx(12)
868 real(kind=
rp) :: curve_data(5,12), x3(3,3,3), y3(3,3,3), z3(3,3,3)
870 real(kind=
rp),
dimension(3),
parameter :: zquad = [-1d0, 0d0,1d0]
871 real(kind=
rp) :: zg(3)
872 real(kind=
rp),
dimension(Xh%lx, Xh%lx, Xh%lx) :: tmp
873 real(kind=
rp) :: jx(xh%lx*3)
874 real(kind=
rp) :: jxt(xh%lx*3), jyt(xh%lx*3), jzt(xh%lx*3)
875 real(kind=
rp) :: w(4*xh%lxyz,2)
876 integer :: j, k, n_edges
877 eindx = [2 , 6 , 8 , 4, &
882 if (msh%gdim .eq. 3)
then
884 call xh3%init(
gll, 3, 3, 3)
887 call xh3%init(
gll, 3, 3)
892 if (curve_type(k) .eq. 4)
then
893 x3(eindx(k),1,1) = curve_data(1,k)
894 y3(eindx(k),1,1) = curve_data(2,k)
895 z3(eindx(k),1,1) = curve_data(3,k)
901 if (msh%gdim .eq. 3)
then
906 call neko_warning(
' m deformation not supported for 2d yet')
914 if (msh%gdim .gt. 2)
then
919 call trsp(jx, xh%lx, jxt, 3)
920 if (msh%gdim .eq. 3)
then
921 call tensr3(tmp, xh%lx, x3, 3, jx, jyt, jzt, w)
922 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
923 call tensr3(tmp, xh%ly, y3, 3, jx, jyt, jzt, w)
924 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
925 call tensr3(tmp, xh%lz, z3, 3, jx, jyt, jzt, w)
926 call copy(z, tmp, xh%lx*xh%ly*xh%lz)
928 call tnsr2d_el(tmp, xh%lx, x3, 3, jx, jyt)
929 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
930 call tnsr2d_el(tmp, xh%ly, y3, 3, jx, jyt)
931 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
943 integer,
intent(in) :: n
944 real(kind=
rp),
intent(inout) :: x(n, n, n)
945 real(kind=
rp),
intent(in) :: zg(n)
946 real(kind=
rp),
intent(inout) :: e(n, n, n)
947 real(kind=
rp),
intent(inout) :: v(n, n, n)
948 integer :: gh_type, ntot, kk, jj, ii, k, j, i
949 real(kind=
xp) :: si, sj, sk, hi, hj, hk
955 do concurrent(i = 1:ntot)
959 do concurrent(i = 1:n, j = 1:n, k = 1:n, &
960 ii = 1:n:n-1, jj = 1:n:n-1, kk = 1:n:n-1)
961 si = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
962 sj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
963 sk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
964 v(i,j,k) = v(i,j,k) + si * sj* sk * x(ii, jj, kk)
967 if (gh_type .eq. 1)
then
968 do concurrent(i = 1:ntot)
976 do concurrent(i = 1:ntot)
982 do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1, kk = 1:n:n-1)
983 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
984 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
985 e(i,j,k) = e(i,j,k) + hj*hk*(x(i, jj, kk) - v(i, jj, kk))
990 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, kk = 1:n:n-1)
991 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
992 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
993 e(i,j,k) = e(i,j,k) + hi*hk*(x(ii, j, kk) - v(ii, j, kk))
998 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, jj = 1:n:n-1)
999 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1000 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1001 e(i,j,k) = e(i,j,k) + hi*hj*(x(ii, jj, k) - v(ii, jj, k))
1004 do concurrent(i = 1:ntot)
1005 e(i,1,1) = e(i,1,1) + v(i,1,1)
1008 if (gh_type .eq. 2)
then
1009 do concurrent(i = 1:ntot)
1017 do concurrent(i = 1:ntot)
1023 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1)
1024 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1025 v(i,j,k) = v(i,j,k) + hi*(x(ii,j,k)-e(ii,j,k))
1031 do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1)
1032 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1033 v(i,j,k) = v(i,j,k) + hj*(x(i, jj, k) - e(i, jj, k))
1039 do concurrent(i = 1:n, j = 1:n, k = 1:n, kk = 1:n:n-1)
1040 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
1041 v(i,j,k) = v(i,j,k) + hk*(x(i, j, kk) - e(i, j, kk))
1044 do concurrent(i = 1:ntot)
1045 v(i,1,1) = v(i,1,1) + e(i,1,1)
1055 integer,
intent(in) :: n
1056 real(kind=
rp),
intent(inout) :: x(n, n)
1057 real(kind=
rp),
intent(in) :: zg(n)
1058 real(kind=
rp),
intent(inout) :: e(n, n)
1059 real(kind=
rp),
intent(inout) :: v(n, n)
1060 integer,
intent(in) :: gh_type
1061 integer :: i,j , jj, ii, ntot
1062 real(kind=
rp) :: si, sj, hi, hj
1072 si = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1073 sj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1074 v(i,j) = v(i,j) + si*sj*x(ii, jj)
1079 if (gh_type .eq. 1)
then
1080 call copy(x, v, ntot)
1092 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1093 e(i,j) = e(i,j) + hj*(x(i, jj) - v(i, jj))
1103 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1104 e(i,j) = e(i,j) + hi*(x(ii,j)-v(ii,j))
1109 call add3(x, e, v, ntot)
1116 integer,
intent(in) :: isid, gdim
1117 type(
space_t),
intent(in) :: Xh
1119 real(kind=
rp),
dimension(5),
intent(in) :: curve_data
1120 real(kind=
rp),
dimension(Xh%lx, Xh%ly, Xh%lz),
intent(inout) :: x, y, z
1121 real(kind=
rp) :: pt1x, pt1y, pt2x, pt2y, pt12x, pt12y
1122 real(kind=
rp) :: radius, dtheta, r, xys
1123 real(kind=
rp) :: theta0, xcenn, ycenn, h(xh%lx, 3, 2)
1124 real(kind=
rp) :: xcrved(xh%lx), ycrved(xh%lx), xs, ys
1125 integer :: isid1, ixt, iyt, izt, ix, itmp
1127 integer(i4),
dimension(6),
parameter :: fcyc_to_sym = [3, 2, 4, 1, 5, 6]
1129 integer(i4),
dimension(12),
parameter :: ecyc_to_sym = [1, 6, 2, 5, 3, 8,&
1130 & 4, 7, 9, 10, 12, 11]
1132 integer,
parameter,
dimension(2, 12) :: edge_nodes = reshape([1, 2, 3, 4, &
1133 & 5, 6, 7, 8, 1, 3, 2, 4, 5, 7, 6, 8, 1, 5, 2, 6, 3, 7, 4, 8], &
1141 itmp = ecyc_to_sym(isid)
1144 pt1x =
element%pts(edge_nodes(1, itmp))%p%x(1)
1145 pt1y =
element%pts(edge_nodes(1, itmp))%p%x(2)
1146 pt2x =
element%pts(edge_nodes(2, itmp))%p%x(1)
1147 pt2y =
element%pts(edge_nodes(2, itmp))%p%x(2)
1149 pt1x =
element%pts(edge_nodes(2, itmp))%p%x(1)
1150 pt1y =
element%pts(edge_nodes(2, itmp))%p%x(2)
1151 pt2x =
element%pts(edge_nodes(1, itmp))%p%x(1)
1152 pt2y =
element%pts(edge_nodes(1, itmp))%p%x(2)
1155 radius = curve_data(1)
1159 xys = sqrt(xs**2 + ys**2)
1161 if (abs(2.0 * radius) <= xys * 1.00001) &
1162 &
call neko_error(
'Radius to small for arced element surface')
1164 dtheta = abs(asin(0.5_xp*xys/radius))
1165 pt12x = (pt1x + pt2x)/2.0
1166 pt12y = (pt1y + pt2y)/2.0
1167 xcenn = pt12x - xs/xys * radius*cos(dtheta)
1168 ycenn = pt12y - ys/xys * radius*cos(dtheta)
1169 theta0 = atan2((pt12y-ycenn), (pt12x-xcenn))
1171 isid1 = mod(isid+4-1, 4)+1
1173 if (radius < 0.0) dtheta = -dtheta
1176 if (isid1 .gt. 2) ixt = xh%lx+1-ix
1178 xcrved(ixt) = xcenn + abs(radius) * cos(theta0 + r*dtheta) &
1179 - ( h(ix,1,1)*pt1x + h(ix,1,2)*pt2x )
1180 ycrved(ixt) = ycenn + abs(radius) * sin(theta0 + r*dtheta) &
1181 - ( h(ix,1,1)*pt1y + h(ix,1,2)*pt2y )
1186 isid1 = fcyc_to_sym(isid1)
1190 if (isid1 .le. 2)
then
1191 call addtnsr(x, h(1, 1, ixt), xcrved, h(1, 3, izt), &
1192 xh%lx, xh%ly, xh%lz)
1193 call addtnsr(y, h(1, 1, ixt), ycrved, h(1, 3, izt), &
1194 xh%lx, xh%ly, xh%lz)
1196 call addtnsr(x, xcrved, h(1, 2, iyt), h(1, 3, izt), &
1197 xh%lx, xh%ly, xh%lz)
1198 call addtnsr(y, ycrved, h(1, 2, iyt), h(1, 3, izt), &
1199 xh%lx, xh%ly, xh%lz)
1204 integer,
intent(in) :: lx, gdim
1205 real(kind=
rp),
intent(inout) :: h(lx, 3, 2)
1206 real(kind=
rp),
intent(in) :: zgml(lx, 3)
1207 integer :: ix, iy, iz
1210 h(ix,1,1) = (1.0_rp - zgml(ix, 1)) * 0.5_rp
1211 h(ix,1,2) = (1.0_rp + zgml(ix, 1)) * 0.5_rp
1215 h(iy,2,1) = (1.0_rp - zgml(iy, 2)) * 0.5_rp
1216 h(iy,2,2) = (1.0_rp + zgml(iy, 2)) * 0.5_rp
1219 if (gdim .eq. 3)
then
1221 h(iz,3,1) = (1.0_rp - zgml(iz, 3)) * 0.5_rp
1222 h(iz,3,2) = (1.0_rp + zgml(iz, 3)) * 0.5_rp
1225 call rone(h(1,3,1), lx)
1226 call rone(h(1,3,2), lx)
Map a Fortran array to a device (allocate and associate)
Copy data between host and device (or device and device)
Device abstraction, common interface for various accelerators.
integer, parameter, public host_to_device
subroutine, public device_free(x_d)
Deallocate memory on the device.
Defines a mapping of the degrees of freedom.
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.
subroutine dofmap_xyzquad(Xh, msh, element, x, y, z, curve_type, curve_data)
subroutine dofmap_generate_xyz(this)
Generate x,y,z-coordinates for all dofs.
subroutine compute_h(h, zgml, gdim, lx)
subroutine arc_surface(isid, curve_data, x, y, z, Xh, element, gdim)
subroutine dofmap_free(this)
Destructor.
subroutine dofmap_init(this, msh, Xh)
Constructor.
subroutine dofmap_number_edges(this)
Assing numbers to dofs on edges.
subroutine dofmap_number_points(this)
Assign numbers to each dofs on points.
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.
pure integer function dofmap_size(this)
Return the total number of dofs in the dofmap, lx*ly*lz*nelv.
subroutine dofmap_number_faces(this)
Assign numbers to dofs on faces.
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,...
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.
Fast diagonalization methods from NEKTON.
subroutine, public fd_weights_full(xi, x, n, m, c)
Compute finite-difference stencil weights for evaluating derivatives up to order at a point.
Defines a hexahedron element.
subroutine, public rone(a, n)
Set all elements to one.
subroutine, public add3(a, b, c, n)
Vector addition .
subroutine, public copy(a, b, n)
Copy a vector .
subroutine, public rzero(a, n)
Zero a real vector.
integer, parameter neko_bcknd_device
integer, parameter, public i8
integer, parameter, public i4
integer, parameter, public xp
integer, parameter, public rp
Global precision used in computations.
Defines a quadrilateral element.
Defines a function space.
integer, parameter, public gll
subroutine, public tensr3(v, nv, u, nu, A, Bt, Ct, w)
Tensor product .
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...
subroutine, public trsp(a, lda, b, ldb)
Transpose of a rectangular tensor .
subroutine, public tnsr2d_el(v, nv, u, nu, A, Bt)
Computes .
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Base type for an element.
The function space for the SEM solution fields.