My Project
Public Member Functions | List of all members
mod_utils::interp_pzonal Interface Reference

Public Member Functions

real(spa) function interp_pzonal_2d_flt (xloc, yloc, i, Field)
 
real(spa) function interp_pzonal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_pzonal_2d_dbl (xloc, yloc, i, Field)
 
real(dp) function interp_pzonal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 

Detailed Description

Definition at line 128 of file mod_utils.f90.

Member Function/Subroutine Documentation

◆ interp_pzonal_2d_dbl()

real(dp) function mod_utils::interp_pzonal::interp_pzonal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), pointer  Field 
)

Definition at line 920 of file mod_utils.f90.

920  USE all_vars, only : a1u,a2u,xc,yc,nbe
921  IMPLICIT NONE
922  REAL(DP) :: FPT
923  INTEGER, INTENT(IN) :: i ! Cell Number
924  REAL(DP), INTENT(IN):: xloc,yloc
925  REAL(DP), POINTER, INTENT(IN) :: FIELD(:)
926 
927  REAL(DP):: X0c, Y0c,F0,Fx,Fy
928  INTEGER :: e1,e2,e3
929  REAL(DP) :: dx_sph
930  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
931 
932 !offset from element center
933 !---------------------------------------------------------------
934  !offset from element center
935  x0c = xloc - xc(i)
936  y0c = yloc - yc(i)
937 !---------------------------------------------------------------
938 !---------------------------------------------------------------
939 
940  !Surrounding Element IDs
941  e1 = nbe(i,1)
942  e2 = nbe(i,2)
943  e3 = nbe(i,3)
944 
945  !interpolate Field to the location
946  fx = a1u(i,1)*field(i)+a1u(i,2)*field(e1)+a1u(i,3)*field(e2)+a1u(i,4)*field(e3)
947  fy = a2u(i,1)*field(i)+a2u(i,2)*field(e1)+a2u(i,3)*field(e2)+a2u(i,4)*field(e3)
948  fpt = field(i) + fx*x0c + fy*y0c
949 
950  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
951 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003

◆ interp_pzonal_2d_flt()

real(spa) function mod_utils::interp_pzonal::interp_pzonal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), pointer  Field 
)

Definition at line 886 of file mod_utils.f90.

886  USE all_vars, only : a1u,a2u,xc,yc,nbe
887  IMPLICIT NONE
888  REAL(SPA) :: FPT
889  INTEGER, INTENT(IN) :: i ! Cell Number
890  REAL(SPA), INTENT(IN):: xloc,yloc
891  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:)
892 
893  REAL(SPA):: X0c, Y0c,F0,Fx,Fy
894  INTEGER :: e1,e2,e3
895  REAL(SPA) :: dx_sph
896  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
897 
898 !offset from element center
899 !---------------------------------------------------------------
900 
901  x0c = xloc - xc(i)
902  y0c = yloc - yc(i)
903 !---------------------------------------------------------------
904 !---------------------------------------------------------------
905  !Surrounding Element IDs
906  e1 = nbe(i,1)
907  e2 = nbe(i,2)
908  e3 = nbe(i,3)
909 
910  !interpolate Field to the location
911  fx = a1u(i,1)*field(i)+a1u(i,2)*field(e1)+a1u(i,3)*field(e2)+a1u(i,4)*field(e3)
912  fy = a2u(i,1)*field(i)+a2u(i,2)*field(e1)+a2u(i,3)*field(e2)+a2u(i,4)*field(e3)
913  fpt = field(i) + fx*x0c + fy*y0c
914 
915  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_2D"
916 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003

◆ interp_pzonal_3d_dbl()

real(dp) function mod_utils::interp_pzonal::interp_pzonal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 1083 of file mod_utils.f90.

1083  USE all_vars, only : a1u,a2u,xc,yc,nbe,kbm2,kbm1,kb,z1,dz1,zz1,dzz1
1084  IMPLICIT NONE
1085  REAL(DP) :: FPT
1086  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of levels(kb,or kbm1)
1087  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
1088  REAL(DP), POINTER, INTENT(IN) :: FIELD(:,:)
1089 
1090  REAL(DP):: X0c, Y0c,F0,Fx,Fy,alpha,F_upper,F_lower
1091  INTEGER :: e1,e2,e3,k1,k2,k
1092 
1093  REAL(DP) :: dx_sph
1094  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1095 
1096 
1097 !offset from element center
1098 !---------------------------------------------------------------
1099  !offset from element center
1100  x0c = xloc - xc(i)
1101  y0c = yloc - yc(i)
1102 !---------------------------------------------------------------
1103 !---------------------------------------------------------------
1104 
1105  !Surrounding Element IDs
1106  e1 = nbe(i,1)
1107  e2 = nbe(i,2)
1108  e3 = nbe(i,3)
1109 
1110 !----determine sigma layers above and below sigloc
1111 
1112  IF(lvls == kbm1) THEN
1113 
1114  k1 = 1
1115  k2 = 1
1116  alpha = -1
1117  if(sigloc < zz1(i,kbm1)) then !!particle near bottom
1118  k1 = kbm1
1119  k2 = kbm1
1120  alpha = -1
1121  else
1122  do k=1,kbm2
1123  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
1124  k1 = k
1125  k2 = k+1
1126  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
1127  exit
1128  endif
1129  end do
1130  end if
1131 
1132  ELSE IF(lvls == kb) THEN
1133 
1134  !surface (default)
1135  k1 = 1
1136  k2 = 1
1137  alpha = -1
1138  !bottom
1139  if(sigloc < z1(i,kb))then
1140  k1 = kb
1141  k2 = kb
1142  alpha = -1
1143  else !intermediate
1144  do k=1,kbm1
1145  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
1146  k1 = k
1147  k2 = k+1
1148  alpha = (z1(i,k)-sigloc)/dz1(i,k)
1149  endif
1150  end do
1151  endif
1152 
1153  ELSE
1154  CALL fatal_error("INTERP_PZONAL_3D: Invalid number of levels passed",&
1155  & "(Must be equal to either KB or KBM1")
1156  END IF
1157 
1158  !interpolate Field to the location
1159  fx = a1u(i,1)*field(i,k1)+a1u(i,2)*field(e1,k1)+a1u(i,3)*field(e2,k1)+a1u(i,4)*field(e3,k1)
1160  fy = a2u(i,1)*field(i,k1)+a2u(i,2)*field(e1,k1)+a2u(i,3)*field(e2,k1)+a2u(i,4)*field(e3,k1)
1161  f_upper = field(i,k1) + fx*x0c + fy*y0c
1162 
1163  IF(k1 == k2) THEN
1164  fpt = f_upper
1165  ELSE
1166 
1167  fx = a1u(i,1)*field(i,k2)+a1u(i,2)*field(e1,k2)+a1u(i,3)*field(e2,k2)+a1u(i,4)*field(e3,k2)
1168  fy = a2u(i,1)*field(i,k2)+a2u(i,2)*field(e1,k2)+a2u(i,3)*field(e2,k2)+a2u(i,4)*field(e3,k2)
1169  f_lower = field(i,k2) + fx*x0c + fy*y0c
1170 
1171  fpt = (alpha)*f_lower + (1.0_dp-alpha)*f_upper
1172  END IF
1173 
1174 
1175  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1176 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer kbm1
Definition: mod_main.f90:65

◆ interp_pzonal_3d_flt()

real(spa) function mod_utils::interp_pzonal::interp_pzonal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 988 of file mod_utils.f90.

988  USE all_vars, only : a1u,a2u,xc,yc,nbe,kbm2,kbm1,kb,z1,dz1,zz1,dzz1
989  IMPLICIT NONE
990  REAL(SPA) :: FPT
991  INTEGER, INTENT(IN) :: i,lvls ! Cell Number, Number of levels(kb,or kbm1)
992  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
993  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:,:)
994 
995  REAL(SPA):: X0c, Y0c,F0,Fx,Fy,alpha,F_upper,F_lower
996  INTEGER :: e1,e2,e3,k1,k2,k
997  REAL(SPA) :: dx_sph
998  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
999 
1000 !offset from element center
1001 !---------------------------------------------------------------
1002  !offset from element center
1003  x0c = xloc - xc(i)
1004  y0c = yloc - yc(i)
1005 !---------------------------------------------------------------
1006 !---------------------------------------------------------------
1007 
1008  !Surrounding Element IDs
1009  e1 = nbe(i,1)
1010  e2 = nbe(i,2)
1011  e3 = nbe(i,3)
1012 
1013 !----determine sigma layers above and below sigloc
1014 
1015  IF(lvls == kbm1) THEN
1016 
1017  k1 = 1
1018  k2 = 1
1019  alpha = -1
1020  if(sigloc < zz1(i,kbm1)) then !!particle near bottom
1021  k1 = kbm1
1022  k2 = kbm1
1023  alpha = -1
1024  else
1025  do k=1,kbm2
1026  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
1027  k1 = k
1028  k2 = k+1
1029  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
1030  exit
1031  endif
1032  end do
1033  end if
1034 
1035  ELSE IF(lvls == kb) THEN
1036 
1037  !surface (default)
1038  k1 = 1
1039  k2 = 1
1040  alpha = -1
1041  !bottom
1042  if(sigloc < z1(i,kb))then
1043  k1 = kb
1044  k2 = kb
1045  alpha = -1
1046  else !intermediate
1047  do k=1,kbm1
1048  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
1049  k1 = k
1050  k2 = k+1
1051  alpha = (z1(i,k)-sigloc)/dz1(i,k)
1052  endif
1053  end do
1054  endif
1055 
1056  ELSE
1057  CALL fatal_error("INTERP_PZONAL_3D: Invalid number of levels passed",&
1058  & "(Must be equal to either KB or KBM1")
1059  END IF
1060 
1061  !interpolate Field to the location
1062  fx = a1u(i,1)*field(i,k1)+a1u(i,2)*field(e1,k1)+a1u(i,3)*field(e2,k1)+a1u(i,4)*field(e3,k1)
1063  fy = a2u(i,1)*field(i,k1)+a2u(i,2)*field(e1,k1)+a2u(i,3)*field(e2,k1)+a2u(i,4)*field(e3,k1)
1064  f_upper = field(i,k1) + fx*x0c + fy*y0c
1065 
1066  IF(k1 == k2) THEN
1067  fpt = f_upper
1068  ELSE
1069 
1070  fx = a1u(i,1)*field(i,k2)+a1u(i,2)*field(e1,k2)+a1u(i,3)*field(e2,k2)+a1u(i,4)*field(e3,k2)
1071  fy = a2u(i,1)*field(i,k2)+a2u(i,2)*field(e1,k2)+a2u(i,3)*field(e2,k2)+a2u(i,4)*field(e3,k2)
1072  f_lower = field(i,k2) + fx*x0c + fy*y0c
1073 
1074  fpt = (alpha)*f_lower + (1.0-alpha)*f_upper
1075  END IF
1076 
1077 
1078  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PZONAL_3D"
1079 
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target a1u
Definition: mod_main.f90:1331
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:,:), allocatable, target a2u
Definition: mod_main.f90:1332
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer kbm1
Definition: mod_main.f90:65

The documentation for this interface was generated from the following file: