236 type(
mesh_t),
pointer :: msh
241 integer(kind=i8) :: num_dofs_edges(3)
242 integer(kind=i8) :: edge_id, edge_offset
243 logical :: shared_dof
249 num_dofs_edges(1) = int(xh%lx - 2,
i8)
250 num_dofs_edges(2) = int(xh%ly - 2,
i8)
251 num_dofs_edges(3) = int(xh%lz - 2,
i8)
252 edge_offset = int(msh%glb_mpts,
i8) + int(1,
i8)
257 select type (ep => msh%elements(i)%e)
262 call ep%edge_id(edge, 1)
263 shared_dof = msh%is_shared(edge)
264 global_id = msh%get_global(edge)
265 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
267 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
268 do concurrent(j = 2:xh%lx - 1)
270 this%dof(k, 1, 1, i) = edge_id + (j-2)
271 this%shared_dof(k, 1, 1, i) = shared_dof
274 do concurrent(j = 2:xh%lx - 1)
276 this%dof(k, 1, 1, i) = edge_id + (j-2)
277 this%shared_dof(k, 1, 1, i) = shared_dof
281 call ep%edge_id(edge, 3)
282 shared_dof = msh%is_shared(edge)
283 global_id = msh%get_global(edge)
284 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
285 if (int(edge%x(1),
i8) .ne. this%dof(1, 1, xh%lz, i))
then
286 do concurrent(j = 2:xh%lx - 1)
288 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
289 this%shared_dof(k, 1, xh%lz, i) = shared_dof
292 do concurrent(j = 2:xh%lx - 1)
294 this%dof(k, 1, xh%lz, i) = edge_id + (j-2)
295 this%shared_dof(k, 1, xh%lz, i) = shared_dof
299 call ep%edge_id(edge, 2)
300 shared_dof = msh%is_shared(edge)
301 global_id = msh%get_global(edge)
302 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
303 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
304 do concurrent(j = 2:xh%lx - 1)
306 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
307 this%shared_dof(k, xh%ly, 1, i) = shared_dof
310 do concurrent(j = 2:xh%lx - 1)
312 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
313 this%shared_dof(k, xh%ly, 1, i) = shared_dof
317 call ep%edge_id(edge, 4)
318 shared_dof = msh%is_shared(edge)
319 global_id = msh%get_global(edge)
320 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
321 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, xh%lz, i))
then
322 do concurrent(j = 2:xh%lx - 1)
324 this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
325 this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
328 do concurrent(j = 2:xh%lx - 1)
330 this%dof(k, xh%ly, xh%lz, i) = edge_id + (j-2)
331 this%shared_dof(k, xh%ly, xh%lz, i) = shared_dof
339 call ep%edge_id(edge, 5)
340 shared_dof = msh%is_shared(edge)
341 global_id = msh%get_global(edge)
342 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
343 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
344 do concurrent(j = 2:xh%ly - 1)
346 this%dof(1, k, 1, i) = edge_id + (j-2)
347 this%shared_dof(1, k, 1, i) = shared_dof
350 do concurrent(j = 2:xh%ly - 1)
352 this%dof(1, k, 1, i) = edge_id + (j-2)
353 this%shared_dof(1, k, 1, i) = shared_dof
357 call ep%edge_id(edge, 7)
358 shared_dof = msh%is_shared(edge)
359 global_id = msh%get_global(edge)
360 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
361 if (int(edge%x(1),
i8) .ne. this%dof(1, 1, xh%lz, i))
then
362 do concurrent(j = 2:xh%ly - 1)
364 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
365 this%shared_dof(1, k, xh%lz, i) = shared_dof
368 do concurrent(j = 2:xh%ly - 1)
370 this%dof(1, k, xh%lz, i) = edge_id + (j-2)
371 this%shared_dof(1, k, xh%lz, i) = shared_dof
375 call ep%edge_id(edge, 6)
376 shared_dof = msh%is_shared(edge)
377 global_id = msh%get_global(edge)
378 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
379 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, 1, 1, i))
then
380 do concurrent(j = 2:xh%ly - 1)
382 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
383 this%shared_dof(xh%lx, k, 1, i) = shared_dof
386 do concurrent(j = 2:xh%ly - 1)
388 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
389 this%shared_dof(xh%lx, k, 1, i) = shared_dof
393 call ep%edge_id(edge, 8)
394 shared_dof = msh%is_shared(edge)
395 global_id = msh%get_global(edge)
396 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
397 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, 1, xh%lz, i))
then
398 do concurrent(j = 2:xh%ly - 1)
400 this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
401 this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
404 do concurrent(j = 2:xh%ly - 1)
406 this%dof(xh%lx, k, xh%lz, i) = edge_id + (j-2)
407 this%shared_dof(xh%lx, k, xh%lz, i) = shared_dof
414 call ep%edge_id(edge, 9)
415 shared_dof = msh%is_shared(edge)
416 global_id = msh%get_global(edge)
417 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
418 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
419 do concurrent(j = 2:xh%lz - 1)
421 this%dof(1, 1, k, i) = edge_id + (j-2)
422 this%shared_dof(1, 1, k, i) = shared_dof
425 do concurrent(j = 2:xh%lz - 1)
427 this%dof(1, 1, k, i) = edge_id + (j-2)
428 this%shared_dof(1, 1, k, i) = shared_dof
432 call ep%edge_id(edge, 10)
433 shared_dof = msh%is_shared(edge)
434 global_id = msh%get_global(edge)
435 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
436 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx,1,1,i))
then
437 do concurrent(j = 2:xh%lz - 1)
439 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
440 this%shared_dof(xh%lx, 1, k, i) = shared_dof
443 do concurrent(j = 2:xh%lz - 1)
445 this%dof(xh%lx, 1, k, i) = edge_id + (j-2)
446 this%shared_dof(xh%lx, 1, k, i) = shared_dof
450 call ep%edge_id(edge, 11)
451 shared_dof = msh%is_shared(edge)
452 global_id = msh%get_global(edge)
453 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
454 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
455 do concurrent(j = 2:xh%lz - 1)
457 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
458 this%shared_dof(1, xh%ly, k, i) = shared_dof
461 do concurrent(j = 2:xh%lz - 1)
463 this%dof(1, xh%ly, k, i) = edge_id + (j-2)
464 this%shared_dof(1, xh%ly, k, i) = shared_dof
468 call ep%edge_id(edge, 12)
469 shared_dof = msh%is_shared(edge)
470 global_id = msh%get_global(edge)
471 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(3)
472 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx, xh%ly, 1, i))
then
473 do concurrent(j = 2:xh%lz - 1)
475 this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
476 this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
479 do concurrent(j = 2:xh%lz - 1)
481 this%dof(xh%lx, xh%ly, k, i) = edge_id + (j-2)
482 this%shared_dof(xh%lx, xh%ly, k, i) = shared_dof
489 call ep%facet_id(edge, 3)
490 shared_dof = msh%is_shared(edge)
491 global_id = msh%get_global(edge)
492 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
494 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
495 do concurrent(j = 2:xh%lx - 1)
497 this%dof(k, 1, 1, i) = edge_id + (j-2)
498 this%shared_dof(k, 1, 1, i) = shared_dof
501 do concurrent(j = 2:xh%lx - 1)
503 this%dof(k, 1, 1, i) = edge_id + (j-2)
504 this%shared_dof(k, 1, 1, i) = shared_dof
508 call ep%facet_id(edge, 4)
509 shared_dof = msh%is_shared(edge)
510 global_id = msh%get_global(edge)
511 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(1)
512 if (int(edge%x(1),
i8) .ne. this%dof(1, xh%ly, 1, i))
then
513 do concurrent(j = 2:xh%lx - 1)
515 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
516 this%shared_dof(k, xh%ly, 1, i) = shared_dof
519 do concurrent(j = 2:xh%lx - 1)
521 this%dof(k, xh%ly, 1, i) = edge_id + (j-2)
522 this%shared_dof(k, xh%ly, 1, i) = shared_dof
529 call ep%facet_id(edge, 1)
530 shared_dof = msh%is_shared(edge)
531 global_id = msh%get_global(edge)
532 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
533 if (int(edge%x(1),
i8) .ne. this%dof(1,1,1,i))
then
534 do concurrent(j = 2:xh%ly - 1)
536 this%dof(1, k, 1, i) = edge_id + (j-2)
537 this%shared_dof(1, k, 1, i) = shared_dof
540 do concurrent(j = 2:xh%ly - 1)
542 this%dof(1, k, 1, i) = edge_id + (j-2)
543 this%shared_dof(1, k, 1, i) = shared_dof
547 call ep%facet_id(edge, 2)
548 shared_dof = msh%is_shared(edge)
549 global_id = msh%get_global(edge)
550 edge_id = edge_offset + int((global_id - 1),
i8) * num_dofs_edges(2)
551 if (int(edge%x(1),
i8) .ne. this%dof(xh%lx,1,1,i))
then
552 do concurrent(j = 2:xh%ly - 1)
554 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
555 this%shared_dof(xh%lx, k, 1, i) = shared_dof
558 do concurrent(j = 2:xh%ly - 1)
560 this%dof(xh%lx, k, 1, i) = edge_id + (j-2)
561 this%shared_dof(xh%lx, k, 1, i) = shared_dof
573 type(
mesh_t),
pointer :: msh
578 integer(kind=i8) :: num_dofs_faces(3)
579 integer(kind=i8) :: facet_offset, facet_id
580 logical :: shared_dof
586 facet_offset = int(msh%glb_mpts,
i8) + &
587 int(msh%glb_meds,
i8) * int(xh%lx-2,
i8) + int(1,
i8)
590 num_dofs_faces(1) = int((xh%ly - 2) * (xh%lz - 2),
i8)
591 num_dofs_faces(2) = int((xh%lx - 2) * (xh%lz - 2),
i8)
592 num_dofs_faces(3) = int((xh%lx - 2) * (xh%ly - 2),
i8)
600 call msh%elements(i)%e%facet_id(face, 1)
601 call msh%elements(i)%e%facet_order(face_order, 1)
602 shared_dof = msh%is_shared(face)
603 global_id = msh%get_global(face)
604 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(1)
605 do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
606 this%dof(1, j, k, i) = &
608 this%shared_dof(1, j, k, i) = shared_dof
611 call msh%elements(i)%e%facet_id(face, 2)
612 call msh%elements(i)%e%facet_order(face_order, 2)
613 shared_dof = msh%is_shared(face)
614 global_id = msh%get_global(face)
615 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(1)
616 do concurrent(j = 2:(xh%ly - 1), k = 2:(xh%lz -1))
617 this%dof(xh%lx, j, k, i) = &
619 this%shared_dof(xh%lx, j, k, i) = shared_dof
626 call msh%elements(i)%e%facet_id(face, 3)
627 call msh%elements(i)%e%facet_order(face_order, 3)
628 shared_dof = msh%is_shared(face)
629 global_id = msh%get_global(face)
630 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(2)
631 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
632 this%dof(j, 1, k, i) = &
634 this%shared_dof(j, 1, k, i) = shared_dof
637 call msh%elements(i)%e%facet_id(face, 4)
638 call msh%elements(i)%e%facet_order(face_order, 4)
639 shared_dof = msh%is_shared(face)
640 global_id = msh%get_global(face)
641 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(2)
642 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%lz - 1))
643 this%dof(j, xh%ly, k, i) = &
645 this%shared_dof(j, xh%ly, k, i) = shared_dof
652 call msh%elements(i)%e%facet_id(face, 5)
653 call msh%elements(i)%e%facet_order(face_order, 5)
654 shared_dof = msh%is_shared(face)
655 global_id = msh%get_global(face)
656 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(3)
657 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
658 this%dof(j, k, 1, i) = &
660 this%shared_dof(j, k, 1, i) = shared_dof
663 call msh%elements(i)%e%facet_id(face, 6)
664 call msh%elements(i)%e%facet_order(face_order, 6)
665 shared_dof = msh%is_shared(face)
666 global_id = msh%get_global(face)
667 facet_id = facet_offset + int((global_id - 1),
i8) * num_dofs_faces(3)
668 do concurrent(j = 2:(xh%lx - 1), k = 2:(xh%ly - 1))
669 this%dof(j, k, xh%lz, i) = &
671 this%shared_dof(j, k, xh%lz, i) = shared_dof
804 type(
mesh_t),
pointer,
intent(in) :: msh
805 type(
space_t),
intent(in) :: Xh
807 real(kind=
rp),
intent(inout) :: x(xh%lx, xh%ly, xh%lz), &
808 y(xh%lx, xh%ly, xh%lz), &
809 z(xh%lx, xh%ly, xh%lz)
810 real(kind=
rp) :: xyzb(2,2,2,3), zgml(xh%lx, 3)
811 real(kind=
rp) :: jx(xh%lx*2)
812 real(kind=
rp) :: jxt(xh%lx*2), jyt(xh%lx*2), jzt(xh%lx*2)
813 real(kind=
rp) :: w(4*xh%lx**3), tmp(xh%lx, xh%lx, xh%lx)
814 real(kind=
rp),
dimension(2),
parameter :: zlin = [-1d0, 1d0]
822 call copy(zgml(1,1), xh%zg(1,1), xh%lx)
823 call copy(zgml(1,2), xh%zg(1,2), xh%ly)
824 if (msh%gdim .gt. 2)
then
825 call copy(zgml(1,3), xh%zg(1,3), xh%lz)
832 if (msh%gdim .gt. 2)
then
837 call trsp(jx, xh%lx, jxt, 2)
839 if (msh%gdim .eq. 2)
then
843 if (msh%gdim .gt. 2)
then
844 do concurrent(j = 1:msh%gdim)
845 xyzb(1,1,1,j) =
element%pts(1)%p%x(j)
846 xyzb(2,1,1,j) =
element%pts(2)%p%x(j)
847 xyzb(1,2,1,j) =
element%pts(3)%p%x(j)
848 xyzb(2,2,1,j) =
element%pts(4)%p%x(j)
850 xyzb(1,1,2,j) =
element%pts(5)%p%x(j)
851 xyzb(2,1,2,j) =
element%pts(6)%p%x(j)
852 xyzb(1,2,2,j) =
element%pts(7)%p%x(j)
853 xyzb(2,2,2,j) =
element%pts(8)%p%x(j)
856 do concurrent(j = 1:msh%gdim)
857 xyzb(1,1,1,j) =
element%pts(1)%p%x(j)
858 xyzb(2,1,1,j) =
element%pts(2)%p%x(j)
859 xyzb(1,2,1,j) =
element%pts(3)%p%x(j)
860 xyzb(2,2,1,j) =
element%pts(4)%p%x(j)
863 if (msh%gdim .eq. 3)
then
864 call tensr3(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt, jzt, w)
865 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
866 call tensr3(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt, jzt, w)
867 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
868 call tensr3(tmp, xh%lz, xyzb(1,1,1,3), 2, jx, jyt, jzt, w)
869 call copy(z, tmp, xh%lx*xh%ly*xh%lz)
871 call tnsr2d_el(tmp, xh%lx, xyzb(1,1,1,1), 2, jx, jyt)
872 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
873 call tnsr2d_el(tmp, xh%ly, xyzb(1,1,1,2), 2, jx, jyt)
874 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
880 type(
mesh_t),
pointer,
intent(in) :: msh
881 type(
space_t),
intent(in) :: Xh
883 real(kind=
rp),
dimension(Xh%lx, Xh%ly, Xh%lz),
intent(inout) :: x, y, z
884 integer :: curve_type(12), eindx(12)
885 real(kind=
rp) :: curve_data(5,12), x3(3,3,3), y3(3,3,3), z3(3,3,3)
887 real(kind=
rp),
dimension(3),
parameter :: zquad = [-1d0, 0d0,1d0]
888 real(kind=
rp) :: zg(3)
889 real(kind=
rp),
dimension(Xh%lx, Xh%lx, Xh%lx) :: tmp
890 real(kind=
rp) :: jx(xh%lx*3)
891 real(kind=
rp) :: jxt(xh%lx*3), jyt(xh%lx*3), jzt(xh%lx*3)
892 real(kind=
rp) :: w(4*xh%lxyz,2)
893 integer :: j, k, n_edges
894 eindx = [2 , 6 , 8 , 4, &
899 if (msh%gdim .eq. 3)
then
901 call xh3%init(
gll, 3, 3, 3)
904 call xh3%init(
gll, 3, 3)
909 if (curve_type(k) .eq. 4)
then
910 x3(eindx(k),1,1) = curve_data(1,k)
911 y3(eindx(k),1,1) = curve_data(2,k)
912 z3(eindx(k),1,1) = curve_data(3,k)
918 if (msh%gdim .eq. 3)
then
923 call neko_warning(
' m deformation not supported for 2d yet')
931 if (msh%gdim .gt. 2)
then
936 call trsp(jx, xh%lx, jxt, 3)
937 if (msh%gdim .eq. 3)
then
938 call tensr3(tmp, xh%lx, x3, 3, jx, jyt, jzt, w)
939 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
940 call tensr3(tmp, xh%ly, y3, 3, jx, jyt, jzt, w)
941 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
942 call tensr3(tmp, xh%lz, z3, 3, jx, jyt, jzt, w)
943 call copy(z, tmp, xh%lx*xh%ly*xh%lz)
945 call tnsr2d_el(tmp, xh%lx, x3, 3, jx, jyt)
946 call copy(x, tmp, xh%lx*xh%ly*xh%lz)
947 call tnsr2d_el(tmp, xh%ly, y3, 3, jx, jyt)
948 call copy(y, tmp, xh%lx*xh%ly*xh%lz)
962 integer,
intent(in) :: n
963 real(kind=
rp),
intent(inout) :: x(n, n, n)
964 real(kind=
rp),
intent(in) :: zg(n)
965 real(kind=
rp),
intent(inout) :: e(n, n, n)
966 real(kind=
rp),
intent(inout) :: v(n, n, n)
967 integer :: gh_type, ntot, kk, jj, ii, k, j, i
968 real(kind=
xp) :: si, sj, sk, hi, hj, hk
974 do concurrent(i = 1:ntot)
978 do concurrent(i = 1:n, j = 1:n, k = 1:n, &
979 ii = 1:n:n-1, jj = 1:n:n-1, kk = 1:n:n-1)
980 si = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
981 sj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
982 sk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
983 v(i,j,k) = v(i,j,k) + si * sj* sk * x(ii, jj, kk)
986 if (gh_type .eq. 1)
then
987 do concurrent(i = 1:ntot)
995 do concurrent(i = 1:ntot)
1001 do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1, kk = 1:n:n-1)
1002 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1003 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
1004 e(i,j,k) = e(i,j,k) + hj*hk*(x(i, jj, kk) - v(i, jj, kk))
1009 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, kk = 1:n:n-1)
1010 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1011 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
1012 e(i,j,k) = e(i,j,k) + hi*hk*(x(ii, j, kk) - v(ii, j, kk))
1017 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1, jj = 1:n:n-1)
1018 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1019 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1020 e(i,j,k) = e(i,j,k) + hi*hj*(x(ii, jj, k) - v(ii, jj, k))
1023 do concurrent(i = 1:ntot)
1024 e(i,1,1) = e(i,1,1) + v(i,1,1)
1027 if (gh_type .eq. 2)
then
1028 do concurrent(i = 1:ntot)
1036 do concurrent(i = 1:ntot)
1042 do concurrent(i = 1:n, j = 1:n, k = 1:n, ii = 1:n:n-1)
1043 hi = 0.5_xp*((n-ii)*(1-zg(i))+(ii-1)*(1+zg(i)))/(n-1)
1044 v(i,j,k) = v(i,j,k) + hi*(x(ii,j,k)-e(ii,j,k))
1050 do concurrent(i = 1:n, j = 1:n, k = 1:n, jj = 1:n:n-1)
1051 hj = 0.5_xp*((n-jj)*(1-zg(j))+(jj-1)*(1+zg(j)))/(n-1)
1052 v(i,j,k) = v(i,j,k) + hj*(x(i, jj, k) - e(i, jj, k))
1058 do concurrent(i = 1:n, j = 1:n, k = 1:n, kk = 1:n:n-1)
1059 hk = 0.5_xp*((n-kk)*(1-zg(k))+(kk-1)*(1+zg(k)))/(n-1)
1060 v(i,j,k) = v(i,j,k) + hk*(x(i, j, kk) - e(i, j, kk))
1063 do concurrent(i = 1:ntot)
1064 v(i,1,1) = v(i,1,1) + e(i,1,1)
1136 integer,
intent(in) :: isid, gdim
1137 type(
space_t),
intent(in) :: Xh
1139 real(kind=
rp),
dimension(5),
intent(in) :: curve_data
1140 real(kind=
rp),
dimension(Xh%lx, Xh%ly, Xh%lz),
intent(inout) :: x, y, z
1141 real(kind=
rp) :: pt1x, pt1y, pt2x, pt2y, pt12x, pt12y
1142 real(kind=
rp) :: radius, dtheta, r, xys
1143 real(kind=
rp) :: theta0, xcenn, ycenn, h(xh%lx, 3, 2)
1144 real(kind=
rp) :: xcrved(xh%lx), ycrved(xh%lx), xs, ys
1145 integer :: isid1, ixt, iyt, izt, ix, itmp
1147 integer(i4),
dimension(6),
parameter :: fcyc_to_sym = [3, 2, 4, 1, 5, 6]
1149 integer(i4),
dimension(12),
parameter :: ecyc_to_sym = [1, 6, 2, 5, 3, 8, &
1150 4, 7, 9, 10, 12, 11]
1152 integer,
parameter,
dimension(2, 12) :: edge_nodes = reshape([1, 2, 3, 4, &
1153 5, 6, 7, 8, 1, 3, 2, 4, 5, 7, 6, 8, 1, 5, 2, 6, 3, 7, 4, 8], &
1161 itmp = ecyc_to_sym(isid)
1164 pt1x =
element%pts(edge_nodes(1, itmp))%p%x(1)
1165 pt1y =
element%pts(edge_nodes(1, itmp))%p%x(2)
1166 pt2x =
element%pts(edge_nodes(2, itmp))%p%x(1)
1167 pt2y =
element%pts(edge_nodes(2, itmp))%p%x(2)
1169 pt1x =
element%pts(edge_nodes(2, itmp))%p%x(1)
1170 pt1y =
element%pts(edge_nodes(2, itmp))%p%x(2)
1171 pt2x =
element%pts(edge_nodes(1, itmp))%p%x(1)
1172 pt2y =
element%pts(edge_nodes(1, itmp))%p%x(2)
1175 radius = curve_data(1)
1179 xys = sqrt(xs**2 + ys**2)
1181 if (abs(2.0 * radius) <= xys * 1.00001) &
1182 &
call neko_error(
'Radius to small for arced element surface')
1184 dtheta = abs(asin(0.5_xp*xys/radius))
1185 pt12x = (pt1x + pt2x)/2.0
1186 pt12y = (pt1y + pt2y)/2.0
1187 xcenn = pt12x - xs/xys * radius*cos(dtheta)
1188 ycenn = pt12y - ys/xys * radius*cos(dtheta)
1189 theta0 = atan2((pt12y-ycenn), (pt12x-xcenn))
1191 isid1 = mod(isid+4-1, 4)+1
1193 if (radius < 0.0) dtheta = -dtheta
1196 if (isid1 .gt. 2) ixt = xh%lx+1-ix
1198 xcrved(ixt) = xcenn + abs(radius) * cos(theta0 + r*dtheta) &
1199 - ( h(ix,1,1)*pt1x + h(ix,1,2)*pt2x )
1200 ycrved(ixt) = ycenn + abs(radius) * sin(theta0 + r*dtheta) &
1201 - ( h(ix,1,1)*pt1y + h(ix,1,2)*pt2y )
1206 isid1 = fcyc_to_sym(isid1)
1210 if (isid1 .le. 2)
then
1211 call addtnsr(x, h(1, 1, ixt), xcrved, h(1, 3, izt), &
1212 xh%lx, xh%ly, xh%lz)
1213 call addtnsr(y, h(1, 1, ixt), ycrved, h(1, 3, izt), &
1214 xh%lx, xh%ly, xh%lz)
1216 call addtnsr(x, xcrved, h(1, 2, iyt), h(1, 3, izt), &
1217 xh%lx, xh%ly, xh%lz)
1218 call addtnsr(y, ycrved, h(1, 2, iyt), h(1, 3, izt), &
1219 xh%lx, xh%ly, xh%lz)