Commit 15adbf52baf8ed530882453ee63cbe9f4d1c2132
1 parent
126a3aed07
Exists in
master
and in
3 other branches
git-svn-id: https://lxsd.femto-st.fr/svn/fvn@24 b657c933-2333-4658-acf2-d3c7c2708721
Showing 2 changed files with 539 additions and 0 deletions Side-by-side Diff
doc/fvn.pdf
No preview for this file type
fvnlib.f90
... | ... | @@ -902,6 +902,545 @@ |
902 | 902 | |
903 | 903 | end subroutine |
904 | 904 | |
905 | + | |
906 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
907 | +! | |
908 | +! Quadratic interpolation of tabulated function of 1,2 or 3 variables | |
909 | +! | |
910 | +!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! | |
911 | + | |
912 | +subroutine fvn_s_find_interval(x,i,xdata,n) | |
913 | + implicit none | |
914 | + ! This routine find the indice i where xdata(i) <= x < xdata(i+1) | |
915 | + ! xdata(n) must contains a set of increasingly ordered values | |
916 | + ! if x < xdata(1) i=0 is returned | |
917 | + ! if x > xdata(n) i=n is returned | |
918 | + ! special case is where x=xdata(n) then n-1 is returned so | |
919 | + ! we will not exclude the upper limit | |
920 | + ! a simple dichotomy method is used | |
921 | + | |
922 | + real(kind=4), intent(in) :: x | |
923 | + real(kind=4), intent(in), dimension(n) :: xdata | |
924 | + integer(kind=4), intent(in) :: n | |
925 | + integer(kind=4), intent(out) :: i | |
926 | + | |
927 | + integer(kind=4) :: imin,imax,imoyen | |
928 | + | |
929 | + ! special case is where x=xdata(n) then n-1 is returned so | |
930 | + ! we will not exclude the upper limit | |
931 | + if (x == xdata(n)) then | |
932 | + i=n-1 | |
933 | + return | |
934 | + end if | |
935 | + | |
936 | + ! if x < xdata(1) i=0 is returned | |
937 | + if (x < xdata(1)) then | |
938 | + i=0 | |
939 | + return | |
940 | + end if | |
941 | + | |
942 | + ! if x > xdata(n) i=n is returned | |
943 | + if (x > xdata(n)) then | |
944 | + i=n | |
945 | + return | |
946 | + end if | |
947 | + | |
948 | + ! here xdata(1) <= x <= xdata(n) | |
949 | + imin=0 | |
950 | + imax=n+1 | |
951 | + | |
952 | + do while((imax-imin) > 1) | |
953 | + imoyen=(imax+imin)/2 | |
954 | + if (x >= xdata(imoyen)) then | |
955 | + imin=imoyen | |
956 | + else | |
957 | + imax=imoyen | |
958 | + end if | |
959 | + end do | |
960 | + | |
961 | + i=imin | |
962 | + | |
963 | +end subroutine | |
964 | + | |
965 | + | |
966 | +subroutine fvn_d_find_interval(x,i,xdata,n) | |
967 | + implicit none | |
968 | + ! This routine find the indice i where xdata(i) <= x < xdata(i+1) | |
969 | + ! xdata(n) must contains a set of increasingly ordered values | |
970 | + ! if x < xdata(1) i=0 is returned | |
971 | + ! if x > xdata(n) i=n is returned | |
972 | + ! special case is where x=xdata(n) then n-1 is returned so | |
973 | + ! we will not exclude the upper limit | |
974 | + ! a simple dichotomy method is used | |
975 | + | |
976 | + real(kind=8), intent(in) :: x | |
977 | + real(kind=8), intent(in), dimension(n) :: xdata | |
978 | + integer(kind=4), intent(in) :: n | |
979 | + integer(kind=4), intent(out) :: i | |
980 | + | |
981 | + integer(kind=4) :: imin,imax,imoyen | |
982 | + | |
983 | + ! special case is where x=xdata(n) then n-1 is returned so | |
984 | + ! we will not exclude the upper limit | |
985 | + if (x == xdata(n)) then | |
986 | + i=n-1 | |
987 | + return | |
988 | + end if | |
989 | + | |
990 | + ! if x < xdata(1) i=0 is returned | |
991 | + if (x < xdata(1)) then | |
992 | + i=0 | |
993 | + return | |
994 | + end if | |
995 | + | |
996 | + ! if x > xdata(n) i=n is returned | |
997 | + if (x > xdata(n)) then | |
998 | + i=n | |
999 | + return | |
1000 | + end if | |
1001 | + | |
1002 | + ! here xdata(1) <= x <= xdata(n) | |
1003 | + imin=0 | |
1004 | + imax=n+1 | |
1005 | + | |
1006 | + do while((imax-imin) > 1) | |
1007 | + imoyen=(imax+imin)/2 | |
1008 | + if (x >= xdata(imoyen)) then | |
1009 | + imin=imoyen | |
1010 | + else | |
1011 | + imax=imoyen | |
1012 | + end if | |
1013 | + end do | |
1014 | + | |
1015 | + i=imin | |
1016 | + | |
1017 | +end subroutine | |
1018 | + | |
1019 | + | |
1020 | +function fvn_s_quad_interpol(x,n,xdata,ydata) | |
1021 | + implicit none | |
1022 | + ! This function evaluate the value of a function defined by a set of points | |
1023 | + ! and values, using a quadratic interpolation | |
1024 | + ! xdata must be increasingly ordered | |
1025 | + ! x must be within xdata(1) and xdata(n) to actually do interpolation | |
1026 | + ! otherwise extrapolation is done | |
1027 | + integer(kind=4), intent(in) :: n | |
1028 | + real(kind=4), intent(in), dimension(n) :: xdata,ydata | |
1029 | + real(kind=4), intent(in) :: x | |
1030 | + real(kind=4) :: fvn_s_quad_interpol | |
1031 | + | |
1032 | + integer(kind=4) :: iinf,base,i,j | |
1033 | + real(kind=4) :: p | |
1034 | + | |
1035 | + call fvn_s_find_interval(x,iinf,xdata,n) | |
1036 | + | |
1037 | + ! Settings for extrapolation | |
1038 | + if (iinf==0) then | |
1039 | + ! TODO -> Lower bound extrapolation warning | |
1040 | + iinf=1 | |
1041 | + end if | |
1042 | + | |
1043 | + if (iinf==n) then | |
1044 | + ! TODO -> Higher bound extrapolation warning | |
1045 | + iinf=n-1 | |
1046 | + end if | |
1047 | + | |
1048 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1049 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1050 | + if (iinf==1) then | |
1051 | + base=0 | |
1052 | + else | |
1053 | + base=iinf-2 | |
1054 | + end if | |
1055 | + | |
1056 | + ! The three points we will use are : | |
1057 | + ! xdata/ydata(base+1),xdata/ydata(base+2),xdata/ydata(base+3) | |
1058 | + | |
1059 | + ! Straight forward Lagrange polynomial | |
1060 | + fvn_s_quad_interpol=0. | |
1061 | + do i=1,3 | |
1062 | + ! polynome i | |
1063 | + p=ydata(base+i) | |
1064 | + do j=1,3 | |
1065 | + if (j /= i) then | |
1066 | + p=p*(x-xdata(base+j))/(xdata(base+i)-xdata(base+j)) | |
1067 | + end if | |
1068 | + end do | |
1069 | + fvn_s_quad_interpol=fvn_s_quad_interpol+p | |
1070 | + end do | |
1071 | + | |
1072 | +end function | |
1073 | + | |
1074 | + | |
1075 | +function fvn_d_quad_interpol(x,n,xdata,ydata) | |
1076 | + implicit none | |
1077 | + ! This function evaluate the value of a function defined by a set of points | |
1078 | + ! and values, using a quadratic interpolation | |
1079 | + ! xdata must be increasingly ordered | |
1080 | + ! x must be within xdata(1) and xdata(n) to actually do interpolation | |
1081 | + ! otherwise extrapolation is done | |
1082 | + integer(kind=4), intent(in) :: n | |
1083 | + real(kind=8), intent(in), dimension(n) :: xdata,ydata | |
1084 | + real(kind=8), intent(in) :: x | |
1085 | + real(kind=8) :: fvn_d_quad_interpol | |
1086 | + | |
1087 | + integer(kind=4) :: iinf,base,i,j | |
1088 | + real(kind=8) :: p | |
1089 | + | |
1090 | + call fvn_d_find_interval(x,iinf,xdata,n) | |
1091 | + | |
1092 | + ! Settings for extrapolation | |
1093 | + if (iinf==0) then | |
1094 | + ! TODO -> Lower bound extrapolation warning | |
1095 | + iinf=1 | |
1096 | + end if | |
1097 | + | |
1098 | + if (iinf==n) then | |
1099 | + ! TODO Higher bound extrapolation warning | |
1100 | + iinf=n-1 | |
1101 | + end if | |
1102 | + | |
1103 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1104 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1105 | + if (iinf==1) then | |
1106 | + base=0 | |
1107 | + else | |
1108 | + base=iinf-2 | |
1109 | + end if | |
1110 | + | |
1111 | + ! The three points we will use are : | |
1112 | + ! xdata/ydata(base+1),xdata/ydata(base+2),xdata/ydata(base+3) | |
1113 | + | |
1114 | + ! Straight forward Lagrange polynomial | |
1115 | + fvn_d_quad_interpol=0. | |
1116 | + do i=1,3 | |
1117 | + ! polynome i | |
1118 | + p=ydata(base+i) | |
1119 | + do j=1,3 | |
1120 | + if (j /= i) then | |
1121 | + p=p*(x-xdata(base+j))/(xdata(base+i)-xdata(base+j)) | |
1122 | + end if | |
1123 | + end do | |
1124 | + fvn_d_quad_interpol=fvn_d_quad_interpol+p | |
1125 | + end do | |
1126 | + | |
1127 | +end function | |
1128 | + | |
1129 | + | |
1130 | +function fvn_s_quad_2d_interpol(x,y,nx,xdata,ny,ydata,zdata) | |
1131 | + implicit none | |
1132 | + ! This function evaluate the value of a two variable function defined by a | |
1133 | + ! set of points and values, using a quadratic interpolation | |
1134 | + ! xdata and ydata must be increasingly ordered | |
1135 | + ! the couple (x,y) must be as x within xdata(1) and xdata(nx) and | |
1136 | + ! y within ydata(1) and ydata(ny) to actually do interpolation | |
1137 | + ! otherwise extrapolation is done | |
1138 | + integer(kind=4), intent(in) :: nx,ny | |
1139 | + real(kind=4), intent(in) :: x,y | |
1140 | + real(kind=4), intent(in), dimension(nx) :: xdata | |
1141 | + real(kind=4), intent(in), dimension(ny) :: ydata | |
1142 | + real(kind=4), intent(in), dimension(nx,ny) :: zdata | |
1143 | + real(kind=4) :: fvn_s_quad_2d_interpol | |
1144 | + | |
1145 | + integer(kind=4) :: ixinf,iyinf,basex,basey,i | |
1146 | + real(kind=4),dimension(3) :: ztmp | |
1147 | + !real(kind=4), external :: fvn_s_quad_interpol | |
1148 | + | |
1149 | + call fvn_s_find_interval(x,ixinf,xdata,nx) | |
1150 | + call fvn_s_find_interval(y,iyinf,ydata,ny) | |
1151 | + | |
1152 | + ! Settings for extrapolation | |
1153 | + if (ixinf==0) then | |
1154 | + ! TODO -> Lower x bound extrapolation warning | |
1155 | + ixinf=1 | |
1156 | + end if | |
1157 | + | |
1158 | + if (ixinf==nx) then | |
1159 | + ! TODO -> Higher x bound extrapolation warning | |
1160 | + ixinf=nx-1 | |
1161 | + end if | |
1162 | + | |
1163 | + if (iyinf==0) then | |
1164 | + ! TODO -> Lower y bound extrapolation warning | |
1165 | + iyinf=1 | |
1166 | + end if | |
1167 | + | |
1168 | + if (iyinf==ny) then | |
1169 | + ! TODO -> Higher y bound extrapolation warning | |
1170 | + iyinf=ny-1 | |
1171 | + end if | |
1172 | + | |
1173 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1174 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1175 | + if (ixinf==1) then | |
1176 | + basex=0 | |
1177 | + else | |
1178 | + basex=ixinf-2 | |
1179 | + end if | |
1180 | + | |
1181 | + if (iyinf==1) then | |
1182 | + basey=0 | |
1183 | + else | |
1184 | + basey=iyinf-2 | |
1185 | + end if | |
1186 | + | |
1187 | + ! First we make 3 interpolations for x at y(base+1),y(base+2),y(base+3) | |
1188 | + ! stored in ztmp(1:3) | |
1189 | + do i=1,3 | |
1190 | + ztmp(i)=fvn_s_quad_interpol(x,nx,xdata,zdata(:,basey+i)) | |
1191 | + end do | |
1192 | + | |
1193 | + ! Then we make an interpolation for y using previous interpolations | |
1194 | + fvn_s_quad_2d_interpol=fvn_s_quad_interpol(y,3,ydata(basey+1:basey+3),ztmp) | |
1195 | +end function | |
1196 | + | |
1197 | + | |
1198 | +function fvn_d_quad_2d_interpol(x,y,nx,xdata,ny,ydata,zdata) | |
1199 | + implicit none | |
1200 | + ! This function evaluate the value of a two variable function defined by a | |
1201 | + ! set of points and values, using a quadratic interpolation | |
1202 | + ! xdata and ydata must be increasingly ordered | |
1203 | + ! the couple (x,y) must be as x within xdata(1) and xdata(nx) and | |
1204 | + ! y within ydata(1) and ydata(ny) to actually do interpolation | |
1205 | + ! otherwise extrapolation is done | |
1206 | + integer(kind=4), intent(in) :: nx,ny | |
1207 | + real(kind=8), intent(in) :: x,y | |
1208 | + real(kind=8), intent(in), dimension(nx) :: xdata | |
1209 | + real(kind=8), intent(in), dimension(ny) :: ydata | |
1210 | + real(kind=8), intent(in), dimension(nx,ny) :: zdata | |
1211 | + real(kind=8) :: fvn_d_quad_2d_interpol | |
1212 | + | |
1213 | + integer(kind=4) :: ixinf,iyinf,basex,basey,i | |
1214 | + real(kind=8),dimension(3) :: ztmp | |
1215 | + !real(kind=8), external :: fvn_d_quad_interpol | |
1216 | + | |
1217 | + call fvn_d_find_interval(x,ixinf,xdata,nx) | |
1218 | + call fvn_d_find_interval(y,iyinf,ydata,ny) | |
1219 | + | |
1220 | + ! Settings for extrapolation | |
1221 | + if (ixinf==0) then | |
1222 | + ! TODO -> Lower x bound extrapolation warning | |
1223 | + ixinf=1 | |
1224 | + end if | |
1225 | + | |
1226 | + if (ixinf==nx) then | |
1227 | + ! TODO -> Higher x bound extrapolation warning | |
1228 | + ixinf=nx-1 | |
1229 | + end if | |
1230 | + | |
1231 | + if (iyinf==0) then | |
1232 | + ! TODO -> Lower y bound extrapolation warning | |
1233 | + iyinf=1 | |
1234 | + end if | |
1235 | + | |
1236 | + if (iyinf==ny) then | |
1237 | + ! TODO -> Higher y bound extrapolation warning | |
1238 | + iyinf=ny-1 | |
1239 | + end if | |
1240 | + | |
1241 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1242 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1243 | + if (ixinf==1) then | |
1244 | + basex=0 | |
1245 | + else | |
1246 | + basex=ixinf-2 | |
1247 | + end if | |
1248 | + | |
1249 | + if (iyinf==1) then | |
1250 | + basey=0 | |
1251 | + else | |
1252 | + basey=iyinf-2 | |
1253 | + end if | |
1254 | + | |
1255 | + ! First we make 3 interpolations for x at y(base+1),y(base+2),y(base+3) | |
1256 | + ! stored in ztmp(1:3) | |
1257 | + do i=1,3 | |
1258 | + ztmp(i)=fvn_d_quad_interpol(x,nx,xdata,zdata(:,basey+i)) | |
1259 | + end do | |
1260 | + | |
1261 | + ! Then we make an interpolation for y using previous interpolations | |
1262 | + fvn_d_quad_2d_interpol=fvn_d_quad_interpol(y,3,ydata(basey+1:basey+3),ztmp) | |
1263 | +end function | |
1264 | + | |
1265 | + | |
1266 | +function fvn_s_quad_3d_interpol(x,y,z,nx,xdata,ny,ydata,nz,zdata,tdata) | |
1267 | + implicit none | |
1268 | + ! This function evaluate the value of a 3 variables function defined by a | |
1269 | + ! set of points and values, using a quadratic interpolation | |
1270 | + ! xdata, ydata and zdata must be increasingly ordered | |
1271 | + ! The triplet (x,y,z) must be within xdata,ydata and zdata to actually | |
1272 | + ! perform an interpolation, otherwise extrapolation is done | |
1273 | + integer(kind=4), intent(in) :: nx,ny,nz | |
1274 | + real(kind=4), intent(in) :: x,y,z | |
1275 | + real(kind=4), intent(in), dimension(nx) :: xdata | |
1276 | + real(kind=4), intent(in), dimension(ny) :: ydata | |
1277 | + real(kind=4), intent(in), dimension(nz) :: zdata | |
1278 | + real(kind=4), intent(in), dimension(nx,ny,nz) :: tdata | |
1279 | + real(kind=4) :: fvn_s_quad_3d_interpol | |
1280 | + | |
1281 | + integer(kind=4) :: ixinf,iyinf,izinf,basex,basey,basez,i,j | |
1282 | + !real(kind=4), external :: fvn_s_quad_interpol,fvn_s_quad_2d_interpol | |
1283 | + real(kind=4),dimension(3,3) :: ttmp | |
1284 | + | |
1285 | + call fvn_s_find_interval(x,ixinf,xdata,nx) | |
1286 | + call fvn_s_find_interval(y,iyinf,ydata,ny) | |
1287 | + call fvn_s_find_interval(z,izinf,zdata,nz) | |
1288 | + | |
1289 | + ! Settings for extrapolation | |
1290 | + if (ixinf==0) then | |
1291 | + ! TODO -> Lower x bound extrapolation warning | |
1292 | + ixinf=1 | |
1293 | + end if | |
1294 | + | |
1295 | + if (ixinf==nx) then | |
1296 | + ! TODO -> Higher x bound extrapolation warning | |
1297 | + ixinf=nx-1 | |
1298 | + end if | |
1299 | + | |
1300 | + if (iyinf==0) then | |
1301 | + ! TODO -> Lower y bound extrapolation warning | |
1302 | + iyinf=1 | |
1303 | + end if | |
1304 | + | |
1305 | + if (iyinf==ny) then | |
1306 | + ! TODO -> Higher y bound extrapolation warning | |
1307 | + iyinf=ny-1 | |
1308 | + end if | |
1309 | + | |
1310 | + if (izinf==0) then | |
1311 | + ! TODO -> Lower z bound extrapolation warning | |
1312 | + izinf=1 | |
1313 | + end if | |
1314 | + | |
1315 | + if (izinf==nz) then | |
1316 | + ! TODO -> Higher z bound extrapolation warning | |
1317 | + izinf=nz-1 | |
1318 | + end if | |
1319 | + | |
1320 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1321 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1322 | + if (ixinf==1) then | |
1323 | + basex=0 | |
1324 | + else | |
1325 | + basex=ixinf-2 | |
1326 | + end if | |
1327 | + | |
1328 | + if (iyinf==1) then | |
1329 | + basey=0 | |
1330 | + else | |
1331 | + basey=iyinf-2 | |
1332 | + end if | |
1333 | + | |
1334 | + if (izinf==1) then | |
1335 | + basez=0 | |
1336 | + else | |
1337 | + basez=izinf-2 | |
1338 | + end if | |
1339 | + | |
1340 | + ! We first make 9 one dimensional interpolation on variable x. | |
1341 | + ! results are stored in ttmp | |
1342 | + do i=1,3 | |
1343 | + do j=1,3 | |
1344 | + ttmp(i,j)=fvn_s_quad_interpol(x,nx,xdata,tdata(:,basey+i,basez+j)) | |
1345 | + end do | |
1346 | + end do | |
1347 | + | |
1348 | + ! We then make a 2 dimensionnal interpolation on variables y and z | |
1349 | + fvn_s_quad_3d_interpol=fvn_s_quad_2d_interpol(y,z, & | |
1350 | + 3,ydata(basey+1:basey+3),3,zdata(basez+1:basez+3),ttmp) | |
1351 | +end function | |
1352 | + | |
1353 | + | |
1354 | +function fvn_d_quad_3d_interpol(x,y,z,nx,xdata,ny,ydata,nz,zdata,tdata) | |
1355 | + implicit none | |
1356 | + ! This function evaluate the value of a 3 variables function defined by a | |
1357 | + ! set of points and values, using a quadratic interpolation | |
1358 | + ! xdata, ydata and zdata must be increasingly ordered | |
1359 | + ! The triplet (x,y,z) must be within xdata,ydata and zdata to actually | |
1360 | + ! perform an interpolation, otherwise extrapolation is done | |
1361 | + integer(kind=4), intent(in) :: nx,ny,nz | |
1362 | + real(kind=8), intent(in) :: x,y,z | |
1363 | + real(kind=8), intent(in), dimension(nx) :: xdata | |
1364 | + real(kind=8), intent(in), dimension(ny) :: ydata | |
1365 | + real(kind=8), intent(in), dimension(nz) :: zdata | |
1366 | + real(kind=8), intent(in), dimension(nx,ny,nz) :: tdata | |
1367 | + real(kind=8) :: fvn_d_quad_3d_interpol | |
1368 | + | |
1369 | + integer(kind=4) :: ixinf,iyinf,izinf,basex,basey,basez,i,j | |
1370 | + !real(kind=8), external :: fvn_d_quad_interpol,fvn_d_quad_2d_interpol | |
1371 | + real(kind=8),dimension(3,3) :: ttmp | |
1372 | + | |
1373 | + call fvn_d_find_interval(x,ixinf,xdata,nx) | |
1374 | + call fvn_d_find_interval(y,iyinf,ydata,ny) | |
1375 | + call fvn_d_find_interval(z,izinf,zdata,nz) | |
1376 | + | |
1377 | + ! Settings for extrapolation | |
1378 | + if (ixinf==0) then | |
1379 | + ! TODO -> Lower x bound extrapolation warning | |
1380 | + ixinf=1 | |
1381 | + end if | |
1382 | + | |
1383 | + if (ixinf==nx) then | |
1384 | + ! TODO -> Higher x bound extrapolation warning | |
1385 | + ixinf=nx-1 | |
1386 | + end if | |
1387 | + | |
1388 | + if (iyinf==0) then | |
1389 | + ! TODO -> Lower y bound extrapolation warning | |
1390 | + iyinf=1 | |
1391 | + end if | |
1392 | + | |
1393 | + if (iyinf==ny) then | |
1394 | + ! TODO -> Higher y bound extrapolation warning | |
1395 | + iyinf=ny-1 | |
1396 | + end if | |
1397 | + | |
1398 | + if (izinf==0) then | |
1399 | + ! TODO -> Lower z bound extrapolation warning | |
1400 | + izinf=1 | |
1401 | + end if | |
1402 | + | |
1403 | + if (izinf==nz) then | |
1404 | + ! TODO -> Higher z bound extrapolation warning | |
1405 | + izinf=nz-1 | |
1406 | + end if | |
1407 | + | |
1408 | + ! The three points we will use are iinf-1,iinf and iinf+1 with the | |
1409 | + ! exception of the first interval, where iinf=1 we will use 1,2 and 3 | |
1410 | + if (ixinf==1) then | |
1411 | + basex=0 | |
1412 | + else | |
1413 | + basex=ixinf-2 | |
1414 | + end if | |
1415 | + | |
1416 | + if (iyinf==1) then | |
1417 | + basey=0 | |
1418 | + else | |
1419 | + basey=iyinf-2 | |
1420 | + end if | |
1421 | + | |
1422 | + if (izinf==1) then | |
1423 | + basez=0 | |
1424 | + else | |
1425 | + basez=izinf-2 | |
1426 | + end if | |
1427 | + | |
1428 | + ! We first make 9 one dimensional interpolation on variable x. | |
1429 | + ! results are stored in ttmp | |
1430 | + do i=1,3 | |
1431 | + do j=1,3 | |
1432 | + ttmp(i,j)=fvn_d_quad_interpol(x,nx,xdata,tdata(:,basey+i,basez+j)) | |
1433 | + end do | |
1434 | + end do | |
1435 | + | |
1436 | + ! We then make a 2 dimensionnal interpolation on variables y and z | |
1437 | + fvn_d_quad_3d_interpol=fvn_d_quad_2d_interpol(y,z, & | |
1438 | + 3,ydata(basey+1:basey+3),3,zdata(basez+1:basez+3),ttmp) | |
1439 | +end function | |
1440 | + | |
1441 | + | |
1442 | + | |
1443 | + | |
905 | 1444 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
906 | 1445 | ! |
907 | 1446 | ! Akima spline interpolation and spline evaluation |