22 #include "sll_assert.h"
23 #include "sll_memory.h"
24 #include "sll_working_precision.h"
25 #include "sll_errors.h"
63 sll_int32 :: num_cells
66 sll_real64 :: delta_eta
82 sll_int32 :: num_cells1
83 sll_int32 :: num_cells2
84 sll_real64 :: eta1_min
85 sll_real64 :: eta1_max
86 sll_real64 :: eta2_min
87 sll_real64 :: eta2_max
88 sll_real64 :: delta_eta1
89 sll_real64 :: delta_eta2
108 sll_int32 :: num_cells1
109 sll_int32 :: num_cells2
110 sll_int32 :: num_cells3
111 sll_real64 :: eta1_min
112 sll_real64 :: eta1_max
113 sll_real64 :: eta2_min
114 sll_real64 :: eta2_max
115 sll_real64 :: eta3_min
116 sll_real64 :: eta3_max
117 sll_real64 :: delta_eta1
118 sll_real64 :: delta_eta2
119 sll_real64 :: delta_eta3
133 sll_int32 :: num_cells1
134 sll_int32 :: num_cells2
135 sll_int32 :: num_cells3
136 sll_int32 :: num_cells4
137 sll_real64 :: eta1_min
138 sll_real64 :: eta1_max
139 sll_real64 :: eta2_min
140 sll_real64 :: eta2_max
141 sll_real64 :: eta3_min
142 sll_real64 :: eta3_max
143 sll_real64 :: eta4_min
144 sll_real64 :: eta4_max
145 sll_real64 :: delta_eta1
146 sll_real64 :: delta_eta2
147 sll_real64 :: delta_eta3
148 sll_real64 :: delta_eta4
162 sll_int32 :: num_cells(6)
163 sll_real64 :: eta_min(6)
164 sll_real64 :: eta_max(6)
165 sll_real64 :: delta_eta(6)
167 sll_real64 :: volume_eta123
168 sll_real64 :: volume_eta456
186 interface operator(*)
189 end interface operator(*)
238 #define TEST_PRESENCE_AND_ASSIGN_VAL( obj, arg, slot, default_val ) \
239 if (
present(arg)) then; \
242 obj%slot = default_val; \
260 sll_int32,
intent(in) :: num_cells
261 sll_real64,
optional,
intent(in) :: eta_min
262 sll_real64,
optional,
intent(in) :: eta_max
265 sll_allocate(m, ierr)
278 sll_int32,
intent(in) :: num_cells
279 sll_real64,
optional,
intent(in) :: eta_min
280 sll_real64,
optional,
intent(in) :: eta_max
282 m%num_cells = num_cells
284 test_presence_and_assign_val(m, eta_min, eta_min, 0.0_f64)
285 test_presence_and_assign_val(m, eta_max, eta_max, 1.0_f64)
287 m%delta_eta = (m%eta_max - m%eta_min)/real(num_cells, f64)
289 if (m%eta_max <= m%eta_min)
then
290 print *,
'ERROR, sll_s_cartesian_mesh_1d_init(): ', &
291 'Problem to construct the mesh 1d '
292 print *,
'because eta_max <= eta_min'
336 sll_real64,
dimension(:),
allocatable :: eta1_node
337 sll_int32 :: num_cells
338 sll_real64 :: eta_min
339 sll_real64 :: delta_eta
343 num_cells = m%num_cells
345 delta_eta = m%delta_eta
346 sll_allocate(eta1_node(num_cells + 1), ierr)
347 do i = 1, num_cells + 1
348 eta1_node(i) = eta_min + real(i - 1, f64)*delta_eta
354 sll_int32,
intent(in) :: i
356 sll_real64 :: eta_min
357 sll_real64 :: delta_eta
359 eta_min = mesh%eta_min
360 delta_eta = mesh%delta_eta
361 res = eta_min + real(i - 1, f64)*delta_eta
366 sll_int32,
intent(in) :: i
368 sll_real64 :: eta_min
369 sll_real64 :: delta_eta
371 eta_min = mesh%eta_min
372 delta_eta = mesh%delta_eta
373 res = eta_min + (real(i - 1, f64) + 0.5_f64)*delta_eta
378 sll_real64,
dimension(:, :),
allocatable :: eta1
379 sll_real64,
dimension(:, :),
allocatable :: eta2
380 sll_int32 :: num_cells1
381 sll_int32 :: num_cells2
382 sll_real64 :: eta1_min
383 sll_real64 :: eta2_min
384 sll_real64 :: delta_eta1
385 sll_real64 :: delta_eta2
390 num_cells1 = m%num_cells1
391 num_cells2 = m%num_cells2
392 eta1_min = m%eta1_min
393 delta_eta1 = m%delta_eta1
394 eta2_min = m%eta2_min
395 delta_eta2 = m%delta_eta2
396 sll_allocate(eta1(num_cells1 + 1, num_cells2 + 1), ierr)
397 sll_allocate(eta2(num_cells1 + 1, num_cells2 + 1), ierr)
398 do j = 1, num_cells2 + 1
399 do i = 1, num_cells1 + 1
400 eta1(i, j) = eta1_min + real(i - 1, f64)*delta_eta1
401 eta2(i, j) = eta2_min + real(j - 1, f64)*delta_eta2
430 sll_int32,
intent(in) :: num_cells1
431 sll_int32,
intent(in) :: num_cells2
432 sll_real64,
optional,
intent(in) :: eta1_min
433 sll_real64,
optional,
intent(in) :: eta1_max
434 sll_real64,
optional,
intent(in) :: eta2_min
435 sll_real64,
optional,
intent(in) :: eta2_max
438 sll_allocate(m, ierr)
473 sll_int32,
intent(in) :: num_cells1
474 sll_int32,
intent(in) :: num_cells2
475 sll_real64,
optional,
intent(in) :: eta1_min
476 sll_real64,
optional,
intent(in) :: eta1_max
477 sll_real64,
optional,
intent(in) :: eta2_min
478 sll_real64,
optional,
intent(in) :: eta2_max
480 m%num_cells1 = num_cells1
481 m%num_cells2 = num_cells2
483 test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
484 test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
485 test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
486 test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
488 m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
489 m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
491 if (m%eta1_max <= m%eta1_min)
then
492 sll_error(
'sll_s_cartesian_mesh_2d_init():',
'Problem to construct the mesh 2d because eta1_max <= eta1_min')
494 if (m%eta2_max <= m%eta2_min)
then
495 sll_error(
'sll_s_cartesian_mesh_2d_init():',
'Problem to construct the mesh 2d because eta2_max <= eta2_min')
501 sll_int32,
intent(in) :: i
502 sll_int32,
intent(in) :: j
504 sll_real64 :: eta1_min
505 sll_real64 :: delta_eta1
511 eta1_min = mesh%eta1_min
512 delta_eta1 = mesh%delta_eta1
513 res = eta1_min + real(i - 1, f64)*delta_eta1
518 sll_int32,
intent(in) :: i
519 sll_int32,
intent(in) :: j
521 sll_real64 :: eta2_min
522 sll_real64 :: delta_eta2
528 eta2_min = mesh%eta2_min
529 delta_eta2 = mesh%delta_eta2
530 res = eta2_min + real(j - 1, f64)*delta_eta2
535 sll_int32,
intent(in) :: i
536 sll_int32,
intent(in) :: j
538 sll_real64 :: eta1_min
539 sll_real64 :: delta_eta1
545 eta1_min = mesh%eta1_min
546 delta_eta1 = mesh%delta_eta1
547 res = eta1_min + (real(i - 1, f64) + 0.5_f64)*delta_eta1
552 sll_int32,
intent(in) :: i
553 sll_int32,
intent(in) :: j
555 sll_real64 :: eta2_min
556 sll_real64 :: delta_eta2
562 eta2_min = mesh%eta2_min
563 delta_eta2 = mesh%delta_eta2
564 res = eta2_min + (real(j - 1, f64) + 0.5_f64)*delta_eta2
569 sll_int32,
intent(in) :: cell_num
571 sll_real64 :: eta1_min
572 sll_real64 :: delta_eta1
575 i = mod(cell_num, mesh%num_cells1)
576 eta1_min = mesh%eta1_min
577 delta_eta1 = mesh%delta_eta1
578 res = eta1_min + (real(i - 1, f64) + 0.5_f64)*delta_eta1
583 sll_int32,
intent(in) :: cell_num
585 sll_real64 :: eta2_min
586 sll_real64 :: delta_eta2
589 j = cell_num/(mesh%num_cells1 + 1) + 1
590 eta2_min = mesh%eta2_min
591 delta_eta2 = mesh%delta_eta2
592 res = eta2_min + (real(j - 1, f64) + 0.5_f64)*delta_eta2
626 sll_int32,
intent(in) :: num_cells1
627 sll_int32,
intent(in) :: num_cells2
628 sll_int32,
intent(in) :: num_cells3
629 sll_real64,
optional,
intent(in) :: eta1_min
630 sll_real64,
optional,
intent(in) :: eta1_max
631 sll_real64,
optional,
intent(in) :: eta2_min
632 sll_real64,
optional,
intent(in) :: eta2_max
633 sll_real64,
optional,
intent(in) :: eta3_min
634 sll_real64,
optional,
intent(in) :: eta3_max
637 sll_allocate(m, ierr)
681 sll_int32,
intent(in) :: num_cells1
682 sll_int32,
intent(in) :: num_cells2
683 sll_int32,
intent(in) :: num_cells3
684 sll_real64,
optional,
intent(in) :: eta1_min
685 sll_real64,
optional,
intent(in) :: eta1_max
686 sll_real64,
optional,
intent(in) :: eta2_min
687 sll_real64,
optional,
intent(in) :: eta2_max
688 sll_real64,
optional,
intent(in) :: eta3_min
689 sll_real64,
optional,
intent(in) :: eta3_max
691 m%num_cells1 = num_cells1
692 m%num_cells2 = num_cells2
693 m%num_cells3 = num_cells3
694 test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
695 test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
696 test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
697 test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
698 test_presence_and_assign_val(m, eta3_min, eta3_min, 0.0_f64)
699 test_presence_and_assign_val(m, eta3_max, eta3_max, 1.0_f64)
700 m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
701 m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
702 m%delta_eta3 = (m%eta3_max - m%eta3_min)/real(num_cells3, f64)
704 if (m%eta1_max <= m%eta1_min)
then
705 print *,
'Problem to construct the mesh 3d '
706 print *,
'because eta1_max <= eta1_min'
708 if (m%eta2_max <= m%eta2_min)
then
709 print *,
'Problem to construct the mesh 3d '
710 print *,
'because eta2_max <= eta2_min'
712 if (m%eta3_max <= m%eta3_min)
then
713 print *,
'Problem to construct the mesh 3d '
714 print *,
'because eta3_max <= eta3_min'
721 sll_int32,
intent(in) :: i1
722 sll_int32,
intent(in) :: i2
723 sll_int32,
intent(in) :: i3
725 sll_real64 :: eta1_min
726 sll_real64 :: delta_eta1
732 eta1_min = mesh%eta1_min
733 delta_eta1 = mesh%delta_eta1
734 res = eta1_min + real(i1 - 1, f64)*delta_eta1
739 sll_int32,
intent(in) :: i1
740 sll_int32,
intent(in) :: i2
741 sll_int32,
intent(in) :: i3
743 sll_real64 :: eta2_min
744 sll_real64 :: delta_eta2
750 eta2_min = mesh%eta2_min
751 delta_eta2 = mesh%delta_eta2
752 res = eta2_min + real(i2 - 1, f64)*delta_eta2
757 sll_int32,
intent(in) :: i1
758 sll_int32,
intent(in) :: i2
759 sll_int32,
intent(in) :: i3
761 sll_real64 :: eta3_min
762 sll_real64 :: delta_eta3
768 eta3_min = mesh%eta3_min
769 delta_eta3 = mesh%delta_eta3
770 res = eta3_min + real(i3 - 1, f64)*delta_eta3
775 sll_int32,
intent(in) :: i1
776 sll_int32,
intent(in) :: i2
777 sll_int32,
intent(in) :: i3
779 sll_real64 :: eta1_min
780 sll_real64 :: delta_eta1
786 eta1_min = mesh%eta1_min
787 delta_eta1 = mesh%delta_eta1
788 res = eta1_min + (real(i1 - 1, f64) + 0.5_f64)*delta_eta1
793 sll_int32,
intent(in) :: i1
794 sll_int32,
intent(in) :: i2
795 sll_int32,
intent(in) :: i3
797 sll_real64 :: eta2_min
798 sll_real64 :: delta_eta2
804 eta2_min = mesh%eta2_min
805 delta_eta2 = mesh%delta_eta2
806 res = eta2_min + (real(i2 - 1, f64) + 0.5_f64)*delta_eta2
811 sll_int32,
intent(in) :: i1
812 sll_int32,
intent(in) :: i2
813 sll_int32,
intent(in) :: i3
815 sll_real64 :: eta3_min
816 sll_real64 :: delta_eta3
822 eta3_min = mesh%eta3_min
823 delta_eta3 = mesh%delta_eta3
824 res = eta3_min + (real(i3 - 1, f64) + 0.5_f64)*delta_eta3
868 sll_int32,
intent(in) :: num_cells1
869 sll_int32,
intent(in) :: num_cells2
870 sll_int32,
intent(in) :: num_cells3
871 sll_int32,
intent(in) :: num_cells4
872 sll_real64,
optional,
intent(in) :: eta1_min
873 sll_real64,
optional,
intent(in) :: eta1_max
874 sll_real64,
optional,
intent(in) :: eta2_min
875 sll_real64,
optional,
intent(in) :: eta2_max
876 sll_real64,
optional,
intent(in) :: eta3_min
877 sll_real64,
optional,
intent(in) :: eta3_max
878 sll_real64,
optional,
intent(in) :: eta4_min
879 sll_real64,
optional,
intent(in) :: eta4_max
882 sll_allocate(m, ierr)
883 m%num_cells1 = num_cells1
884 m%num_cells2 = num_cells2
885 m%num_cells3 = num_cells3
886 m%num_cells4 = num_cells4
887 test_presence_and_assign_val(m, eta1_min, eta1_min, 0.0_f64)
888 test_presence_and_assign_val(m, eta1_max, eta1_max, 1.0_f64)
889 test_presence_and_assign_val(m, eta2_min, eta2_min, 0.0_f64)
890 test_presence_and_assign_val(m, eta2_max, eta2_max, 1.0_f64)
891 test_presence_and_assign_val(m, eta3_min, eta3_min, 0.0_f64)
892 test_presence_and_assign_val(m, eta3_max, eta3_max, 1.0_f64)
893 test_presence_and_assign_val(m, eta4_min, eta4_min, 0.0_f64)
894 test_presence_and_assign_val(m, eta4_max, eta4_max, 1.0_f64)
895 m%delta_eta1 = (m%eta1_max - m%eta1_min)/real(num_cells1, f64)
896 m%delta_eta2 = (m%eta2_max - m%eta2_min)/real(num_cells2, f64)
897 m%delta_eta3 = (m%eta3_max - m%eta3_min)/real(num_cells3, f64)
898 m%delta_eta4 = (m%eta4_max - m%eta4_min)/real(num_cells4, f64)
900 if (m%eta1_max <= m%eta1_min)
then
901 print *,
'Problem to construct the mesh 4d '
902 print *,
'because eta1_max <= eta1_min'
904 if (m%eta2_max <= m%eta2_min)
then
905 print *,
'Problem to construct the mesh 4d '
906 print *,
'because eta2_max <= eta2_min'
908 if (m%eta3_max <= m%eta3_min)
then
909 print *,
'Problem to construct the mesh 4d '
910 print *,
'because eta3_max <= eta3_min'
912 if (m%eta4_max <= m%eta4_min)
then
913 print *,
'Problem to construct the mesh 4d '
914 print *,
'because eta4_max <= eta4_min'
920 sll_int32,
intent(in) :: i1
921 sll_int32,
intent(in) :: i2
922 sll_int32,
intent(in) :: i3
923 sll_int32,
intent(in) :: i4
925 sll_real64 :: eta1_min
926 sll_real64 :: delta_eta1
929 dummy = i1 + i2 + i3 + i4
932 eta1_min = mesh%eta1_min
933 delta_eta1 = mesh%delta_eta1
934 res = eta1_min + real(i1 - 1, f64)*delta_eta1
939 sll_int32,
intent(in) :: i1
940 sll_int32,
intent(in) :: i2
941 sll_int32,
intent(in) :: i3
942 sll_int32,
intent(in) :: i4
944 sll_real64 :: eta2_min
945 sll_real64 :: delta_eta2
948 dummy = i1 + i2 + i3 + i4
951 eta2_min = mesh%eta2_min
952 delta_eta2 = mesh%delta_eta2
953 res = eta2_min + real(i2 - 1, f64)*delta_eta2
958 sll_int32,
intent(in) :: i1
959 sll_int32,
intent(in) :: i2
960 sll_int32,
intent(in) :: i3
961 sll_int32,
intent(in) :: i4
963 sll_real64 :: eta3_min
964 sll_real64 :: delta_eta3
967 dummy = i1 + i2 + i3 + i4
970 eta3_min = mesh%eta3_min
971 delta_eta3 = mesh%delta_eta3
972 res = eta3_min + real(i3 - 1, f64)*delta_eta3
977 sll_int32,
intent(in) :: i1
978 sll_int32,
intent(in) :: i2
979 sll_int32,
intent(in) :: i3
980 sll_int32,
intent(in) :: i4
982 sll_real64 :: eta4_min
983 sll_real64 :: delta_eta4
986 dummy = i1 + i2 + i3 + i4
989 eta4_min = mesh%eta4_min
990 delta_eta4 = mesh%delta_eta4
991 res = eta4_min + real(i4 - 1, f64)*delta_eta4
996 sll_int32,
intent(in) :: i1
997 sll_int32,
intent(in) :: i2
998 sll_int32,
intent(in) :: i3
999 sll_int32,
intent(in) :: i4
1001 sll_real64 :: eta1_min
1002 sll_real64 :: delta_eta1
1005 dummy = i1 + i2 + i3 + i4
1008 eta1_min = mesh%eta1_min
1009 delta_eta1 = mesh%delta_eta1
1010 res = eta1_min + (real(i1 - 1, f64) + 0.5_f64)*delta_eta1
1015 sll_int32,
intent(in) :: i1
1016 sll_int32,
intent(in) :: i2
1017 sll_int32,
intent(in) :: i3
1018 sll_int32,
intent(in) :: i4
1020 sll_real64 :: eta2_min
1021 sll_real64 :: delta_eta2
1024 dummy = i1 + i2 + i3 + i4
1027 eta2_min = mesh%eta2_min
1028 delta_eta2 = mesh%delta_eta2
1029 res = eta2_min + (real(i2 - 1, f64) + 0.5_f64)*delta_eta2
1034 sll_int32,
intent(in) :: i1
1035 sll_int32,
intent(in) :: i2
1036 sll_int32,
intent(in) :: i3
1037 sll_int32,
intent(in) :: i4
1039 sll_real64 :: eta3_min
1040 sll_real64 :: delta_eta3
1043 dummy = i1 + i2 + i3 + i4
1046 eta3_min = mesh%eta3_min
1047 delta_eta3 = mesh%delta_eta3
1048 res = eta3_min + (real(i3 - 1, f64) + 0.5_f64)*delta_eta3
1053 sll_int32,
intent(in) :: i1
1054 sll_int32,
intent(in) :: i2
1055 sll_int32,
intent(in) :: i3
1056 sll_int32,
intent(in) :: i4
1058 sll_real64 :: eta4_min
1059 sll_real64 :: delta_eta4
1062 dummy = i1 + i2 + i3 + i4
1065 eta4_min = mesh%eta4_min
1066 delta_eta4 = mesh%delta_eta4
1067 res = eta4_min + (real(i4 - 1, f64) + 0.5_f64)*delta_eta4
1076 write (*,
"(/,(a))")
'1D mesh : num_cell eta_min eta_max delta_eta'
1077 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells, &
1090 write (*,
"(/,(a))")
'2D mesh : num_cell eta_min eta_max delta_eta'
1091 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1095 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1107 write (*,
"(/,(a))")
'3D mesh : num_cell eta_min eta_max delta_eta'
1108 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1112 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1116 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells3, &
1128 write (*,
"(/,(a))")
'4D mesh : num_cell eta_min eta_max delta_eta'
1129 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells1, &
1133 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells2, &
1137 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells3, &
1141 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells4, &
1158 sll_assert(mesh%num_cells > 0)
1172 sll_assert(mesh%num_cells1 > 0)
1186 sll_assert(mesh%num_cells1 > 0)
1200 sll_assert(mesh%num_cells1 > 0)
1207 sll_real64,
dimension(mesh%num_cells + 1) ::nodes
1210 nknots = mesh%num_cells + 1
1212 nodes(idx) = mesh%eta_min + real(idx - 1, f64)*mesh%delta_eta
1214 nodes(1) = mesh%eta_min
1215 nodes(nknots) = mesh%eta_max
1224 sll_real64,
dimension(:),
intent(in) ::point
1225 sll_int32,
dimension(size(point)) :: cell
1226 cell = floor((point - mesh%eta_min)/mesh%delta_eta) + 1
1228 where (cell == mesh%num_cells + 1) cell = mesh%num_cells
1238 sll_int32,
intent(in) :: cell
1239 sll_real64,
dimension(2) :: margin
1241 margin(1) = mesh%eta_min + real(cell - 1, f64)*mesh%delta_eta
1242 margin(2) = mesh%eta_min + real(cell, f64)*mesh%delta_eta
1250 sll_int32 :: num_nodes
1251 num_nodes = mesh%num_cells + 1
1258 sll_real64 :: length
1259 length = mesh%eta_max - mesh%eta_min
1260 sll_assert(length >= 0)
1268 area = (mesh%eta1_max - mesh%eta1_min)* &
1269 (mesh%eta2_max - mesh%eta2_min)
1271 sll_assert(area >= 0)
1279 area = (mesh%eta1_max - mesh%eta1_min)* &
1280 (mesh%eta2_max - mesh%eta2_min)* &
1281 (mesh%eta3_max - mesh%eta3_min)
1283 sll_assert(area >= 0)
1291 area = (mesh%eta1_max - mesh%eta1_min)* &
1292 (mesh%eta2_max - mesh%eta2_min)* &
1293 (mesh%eta3_max - mesh%eta3_min)* &
1294 (mesh%eta4_max - mesh%eta4_min)
1296 sll_assert(area >= 0)
1304 sll_real64,
intent(in) :: point
1305 sll_real64 :: posmod
1307 posmod = modulo(point - mesh%eta_min, mesh%eta_max - mesh%eta_min) + mesh%eta_min
1309 sll_assert(posmod >= mesh%eta_min)
1310 sll_assert(posmod <= mesh%eta_max)
1326 sll_int32,
intent(in) :: num_cells(6)
1327 sll_real64,
intent(in) :: eta_min(6)
1328 sll_real64,
intent(in) :: eta_max(6)
1330 this%num_cells = num_cells
1331 this%eta_min = eta_min
1332 this%eta_max = eta_max
1334 this%delta_eta = (eta_max - eta_min)/real(num_cells, f64)
1336 this%volume_eta123 = product(this%delta_eta(1:3))
1337 this%volume_eta456 = product(this%delta_eta(4:6))
1338 this%volume = this%volume_eta123*this%volume_eta456
1348 write (*,
"(/,(a))")
'6D mesh : num_cell eta_min eta_max delta_eta'
1349 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(1), &
1353 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(2), &
1357 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(3), &
1361 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(4), &
1365 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(5), &
1369 write (*,
"(10x,(i4,1x),3(g13.3,1x))") mesh%num_cells(6), &
1382 #undef TEST_PRESENCE_AND_ASSIGN_VAL
Deallocates memory for the cartesian mesh.
Get node positions array.
allocates the memory space for a new cartesian mesh on the heap, initializes it with the given argume...
Cartesian mesh basic types.
real(kind=f64) function, dimension(2) cell_margin_cartesian_mesh_1d(mesh, cell)
Returns the margin (a,b) for a given cell.
real(kind=f64) function eta3_node_3d(mesh, i1, i2, i3)
real(kind=f64) function eta1_node_1d(mesh, i)
subroutine display_cartesian_mesh_1d(mesh)
display contents of a 1D cartesian mesh. Recommended access through the generic interface sll_o_displ...
integer(kind=i32) function, dimension(size(point)) cell_cartesian_mesh_1d(mesh, point)
Returns cell number(s) for given point(s) in cartesian mesh.
subroutine sll_s_cartesian_mesh_4d_free(mesh)
deallocates memory for the 4D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function eta1_cell_2d_one_arg(mesh, cell_num)
real(kind=f64) function eta1_cell_4d(mesh, i1, i2, i3, i4)
subroutine sll_s_cartesian_mesh_2d_free(mesh)
deallocates memory for the 2D cartesian mesh. Recommended access through the generic interface delete...
subroutine display_cartesian_mesh_6d(mesh)
display contents of a 6d cartesian mesh. Recommended access through the generic interface sll_display...
subroutine display_cartesian_mesh_4d(mesh)
display contents of a 4d cartesian mesh. Recommended access through the generic interface sll_o_displ...
subroutine, public sll_s_cartesian_mesh_3d_free(mesh)
deallocates memory for the 3D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function eta2_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta2_node_2d(mesh, i, j)
subroutine, public sll_s_cartesian_mesh_3d_init(m, num_cells1, num_cells2, num_cells3, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max)
allocates the memory space for a new 3D cartesian mesh, initializes it with the given arguments and r...
subroutine sll_s_cartesian_mesh_6d_free(this)
destructor of 6D cartesian mesh
real(kind=f64) function, dimension(mesh%num_cells+1) nodes_cartesian_mesh_1d(mesh)
Returns all nodes for the 1D cartesian mesh.
type(sll_t_cartesian_mesh_2d) function, pointer, public sll_f_new_cartesian_mesh_2d(num_cells1, num_cells2, eta1_min, eta1_max, eta2_min, eta2_max)
allocates the memory space for a new 2D cartesian mesh on the heap, initializes it with the given arg...
real(kind=f64) function eta2_cell_2d_two_arg(mesh, i, j)
real(kind=f64) function eta3_cell_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta4_node_4d(mesh, i1, i2, i3, i4)
subroutine get_node_positions_2d(m, eta1, eta2)
subroutine sll_s_cartesian_mesh_1d_free(mesh)
deallocates memory for the 1D cartesian mesh. Recommended access through the generic interface delete...
real(kind=f64) function period_cartesian_mesh_1d(mesh, point)
Returns coordinate in a periodic mesh.
type(sll_t_cartesian_mesh_4d) function, pointer, public sll_f_new_cartesian_mesh_4d(num_cells1, num_cells2, num_cells3, num_cells4, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max, eta4_min, eta4_max)
allocates the memory space for a new 4D cartesian mesh on the heap,
type(sll_t_cartesian_mesh_3d) function, pointer, public sll_f_new_cartesian_mesh_3d(num_cells1, num_cells2, num_cells3, eta1_min, eta1_max, eta2_min, eta2_max, eta3_min, eta3_max)
allocates the memory space for a new 3D cartesian mesh on the heap, initializes it with the given arg...
subroutine sll_s_cartesian_mesh_6d_init(this, num_cells, eta_min, eta_max)
initialize a 6D cartesian mesh
real(kind=f64) function area_cartesian_mesh_4d(mesh)
Returns the area size for a 3D cartesian mesh.
real(kind=f64) function eta1_cell_2d_two_arg(mesh, i, j)
real(kind=f64) function length_cartesian_mesh_1d(mesh)
Returns the interval length for a 1D cartesian mesh.
integer(kind=i32) function num_nodes_cartesian_mesh_1d(mesh)
Returns the number of nodes for the 1D cartesian mesh.
type(sll_t_cartesian_mesh_1d) function, pointer, public sll_f_new_cartesian_mesh_1d(num_cells, eta_min, eta_max)
allocates the memory space for a new 1D cartesian mesh on the heap, initializes it with the given arg...
real(kind=f64) function eta2_node_3d(mesh, i1, i2, i3)
subroutine display_cartesian_mesh_3d(mesh)
display contents of a 3d cartesian mesh. Recommended access through the generic interface sll_o_displ...
subroutine, public sll_s_cartesian_mesh_1d_init(m, num_cells, eta_min, eta_max)
Initializes a previously allocated 1D cartesian mesh object.
real(kind=f64) function eta1_cell_1d(mesh, i)
real(kind=f64) function eta3_node_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta2_cell_2d_one_arg(mesh, cell_num)
real(kind=f64) function eta1_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta1_node_4d(mesh, i1, i2, i3, i4)
subroutine, public sll_s_cartesian_mesh_2d_init(m, num_cells1, num_cells2, eta1_min, eta1_max, eta2_min, eta2_max)
initializes a cartesian mesh 2D object that has been already allocated.
real(kind=f64) function area_cartesian_mesh_2d(mesh)
Returns the area size for a 2D cartesian mesh.
type(sll_t_cartesian_mesh_2d) function, pointer, public sll_f_tensor_product_1d_1d(m_a, m_b)
Create a 2d mesh from two 1d meshes.
real(kind=f64) function eta2_node_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function area_cartesian_mesh_3d(mesh)
Returns the area size for a 3D cartesian mesh.
subroutine display_cartesian_mesh_2d(mesh)
display contents of a 2d cartesian mesh. Recommended access through the generic interface sll_o_displ...
real(kind=f64) function eta1_node_2d(mesh, i, j)
real(kind=f64) function eta4_cell_4d(mesh, i1, i2, i3, i4)
real(kind=f64) function eta1_node_3d(mesh, i1, i2, i3)
type(sll_t_cartesian_mesh_4d) function, pointer tensor_product_2d_2d(m_a, m_b)
Create a 4d mesh from two 2d meshes.
real(kind=f64) function eta3_cell_3d(mesh, i1, i2, i3)
real(kind=f64) function eta2_cell_4d(mesh, i1, i2, i3, i4)
subroutine get_node_positions_1d(m, eta1_node)
2d cartesian mesh pointer