286 type(
mesh_t),
pointer :: msh
291 integer(kind=i8) :: num_dofs_edges(3)
292 integer(kind=i8) :: edge_id, edge_offset
293 logical :: shared_dof
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)
307 select type (ep => msh%elements(i)%e)
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)
317 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
318 do concurrent(j = 2:xh%lx - 1)
320 this%dof(k, 1, 1, i) = edge_id + (j-2)
321 this%shared_dof(k, 1, 1, i) = shared_dof
324 do concurrent(j = 2:xh%lx - 1)
326 this%dof(k, 1, 1, i) = edge_id + (j-2)
327 this%shared_dof(k, 1, 1, i) = shared_dof
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)
338 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
339 this%shared_dof(k, 1, xh%lz, i) = shared_dof
342 do concurrent(j = 2:xh%lx - 1)
344 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
345 this%shared_dof(k, 1, xh%lz, i) = shared_dof
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)
356 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
357 this%shared_dof(k, xh%ly, 1, i) = shared_dof
360 do concurrent(j = 2:xh%lx - 1)
362 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
363 this%shared_dof(k, xh%ly, 1, i) = shared_dof
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)
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
378 do concurrent(j = 2:xh%lx - 1)
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
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)
396 this%dof(1, k, 1, i) = edge_id + (j-2)
397 this%shared_dof(1, k, 1, i) = shared_dof
400 do concurrent(j = 2:xh%ly - 1)
402 this%dof(1, k, 1, i) = edge_id + (j-2)
403 this%shared_dof(1, k, 1, i) = shared_dof
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)
414 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
415 this%shared_dof(1, k, xh%lz, i) = shared_dof
418 do concurrent(j = 2:xh%ly - 1)
420 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
421 this%shared_dof(1, k, xh%lz, i) = shared_dof
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)
432 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
433 this%shared_dof(xh%lx, k, 1, i) = shared_dof
436 do concurrent(j = 2:xh%ly - 1)
438 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
439 this%shared_dof(xh%lx, k, 1, i) = shared_dof
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)
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
454 do concurrent(j = 2:xh%ly - 1)
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
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)
471 this%dof(1, 1, k, i) = edge_id + (j-2)
472 this%shared_dof(1, 1, k, i) = shared_dof
475 do concurrent(j = 2:xh%lz - 1)
477 this%dof(1, 1, k, i) = edge_id + (j-2)
478 this%shared_dof(1, 1, k, i) = shared_dof
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)
489 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
490 this%shared_dof(xh%lx, 1, k, i) = shared_dof
493 do concurrent(j = 2:xh%lz - 1)
495 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
496 this%shared_dof(xh%lx, 1, k, i) = shared_dof
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)
507 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
508 this%shared_dof(1, xh%ly, k, i) = shared_dof
511 do concurrent(j = 2:xh%lz - 1)
513 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
514 this%shared_dof(1, xh%ly, k, i) = shared_dof
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)
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
529 do concurrent(j = 2:xh%lz - 1)
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
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)
544 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
545 do concurrent(j = 2:xh%lx - 1)
547 this%dof(k, 1, 1, i) = edge_id + (j-2)
548 this%shared_dof(k, 1, 1, i) = shared_dof
551 do concurrent(j = 2:xh%lx - 1)
553 this%dof(k, 1, 1, i) = edge_id + (j-2)
554 this%shared_dof(k, 1, 1, i) = shared_dof
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)
565 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
566 this%shared_dof(k, xh%ly, 1, i) = shared_dof
569 do concurrent(j = 2:xh%lx - 1)
571 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
572 this%shared_dof(k, xh%ly, 1, i) = shared_dof
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)
586 this%dof(1, k, 1, i) = edge_id + (j-2)
587 this%shared_dof(1, k, 1, i) = shared_dof
590 do concurrent(j = 2:xh%ly - 1)
592 this%dof(1, k, 1, i) = edge_id + (j-2)
593 this%shared_dof(1, k, 1, i) = shared_dof
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)
604 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
605 this%shared_dof(xh%lx, k, 1, i) = shared_dof
608 do concurrent(j = 2:xh%ly - 1)
610 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
611 this%shared_dof(xh%lx, k, 1, i) = shared_dof
623 type(
mesh_t),
pointer :: msh
628 integer(kind=i8) :: num_dofs_faces(3)
629 integer(kind=i8) :: facet_offset, facet_id
630 logical :: shared_dof
636 facet_offset = int(msh%glb_mpts,
i8) + &
637 int(msh%glb_meds,
i8) * int(xh%lx-2,
i8) + int(1,
i8)
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)
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)
658 facet_id, j, k, xh%lz, xh%ly)
659 this%shared_dof(1, j, k, i) = shared_dof
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)
671 facet_id, j, k, xh%lz, xh%ly)
672 this%shared_dof(xh%lx, j, k, i) = shared_dof
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)
688 facet_id, k, j, xh%lz, xh%lx)
689 this%shared_dof(j, 1, k, i) = shared_dof
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)
701 facet_id, k, j, xh%lz, xh%lx)
702 this%shared_dof(j, xh%ly, k, i) = shared_dof
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)
718 facet_id, k, j, xh%ly, xh%lx)
719 this%shared_dof(j, k, 1, i) = shared_dof
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)
731 facet_id, k, j, xh%lz, xh%lx)
732 this%shared_dof(j, k, xh%lz, i) = shared_dof
866 type(
mesh_t),
pointer,
intent(in) :: msh
867 type(
space_t),
intent(in) :: Xh
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]
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)
894 if (msh%gdim .gt. 2)
then
899 call trsp(jx, xh%lx, jxt, 2)
901 if (msh%gdim .eq. 2)
then
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)
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)
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)
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)
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)
942 type(
mesh_t),
pointer,
intent(in) :: msh
943 type(
space_t),
intent(in) :: Xh
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)
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, &
961 if (msh%gdim .eq. 3)
then
963 call xh3%init(
gll, 3, 3, 3)
966 call xh3%init(
gll, 3, 3)
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)
980 if (msh%gdim .eq. 3)
then
985 call neko_warning(
' m deformation not supported for 2d yet')
993 if (msh%gdim .gt. 2)
then
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)
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)
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
1036 do concurrent(i = 1:ntot)
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)
1048 if (gh_type .eq. 1)
then
1049 do concurrent(i = 1:ntot)
1057 do concurrent(i = 1:ntot)
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))
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))
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))
1085 do concurrent(i = 1:ntot)
1086 e(i,1,1) = e(i,1,1) + v(i,1,1)
1089 if (gh_type .eq. 2)
then
1090 do concurrent(i = 1:ntot)
1098 do concurrent(i = 1:ntot)
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))
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))
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))
1125 do concurrent(i = 1:ntot)
1126 v(i,1,1) = v(i,1,1) + e(i,1,1)
1198 integer,
intent(in) :: isid, gdim
1199 type(
space_t),
intent(in) :: Xh
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
1209 integer(i4),
dimension(6),
parameter :: fcyc_to_sym = [3, 2, 4, 1, 5, 6]
1211 integer(i4),
dimension(12),
parameter :: ecyc_to_sym = [1, 6, 2, 5, 3, 8, &
1212 4, 7, 9, 10, 12, 11]
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], &
1223 itmp = ecyc_to_sym(isid)
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)
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)
1237 radius = curve_data(1)
1241 xys = sqrt(xs**2 + ys**2)
1243 if (abs(2.0 * radius) <= xys * 1.00001) &
1244 &
call neko_error(
'Radius to small for arced element surface')
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))
1253 isid1 = mod(isid+4-1, 4)+1
1255 if (radius < 0.0) dtheta = -dtheta
1258 if (isid1 .gt. 2) ixt = xh%lx+1-ix
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 )
1268 isid1 = fcyc_to_sym(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)
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)