FE-Project
Loading...
Searching...
No Matches
mod_atmos_vars Module Reference

module Atmosphere / Physics / Variables More...

Data Types

type  atmosvars
 Derived type to manage variables with atmospheric component. More...

Functions/Subroutines

subroutine, public atmosvars_getlocalmeshprgvar (domid, mesh, prgvars_list, auxvars_list, varid, var, dens_hyd, pres_hyd, lcmesh3d)
subroutine, public atmosvars_getlocalmeshprgvars (domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
subroutine, public atmosvars_getlocalmeshsfcvar (domid, mesh, auxvars2d_list, prec, prec_engi, lcmesh2d)
subroutine, public atmosvars_getlocalmeshqtrcvar (domid, mesh, trcvars_list, varid, var, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrc_qv (domid, mesh, trcvars_list, forcing_list, var, var_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcvarlist (domid, mesh, trcvars_list, varid_s, var_list, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphyauxvars (domid, mesh, phyauxvars_list, pres, pt, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphytends (domid, mesh, phytends_list, dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, rhoh_p, rhoq_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcphytend (domid, mesh, phytends_list, qtrcid, rhoq_tp)

Variables

real(rp), dimension(prgvar_num), parameter progvars_check_min = (/ -1.0_RP, -100.0_RP, -200.0_RP, -200.0_RP, -200.0_RP /)
integer, parameter, public atmos_auxvars2d_prec_id = 1
integer, parameter, public atmos_auxvars2d_prec_engi_id = 2
integer, parameter, public atmos_auxvars2d_num = 2
type(variableinfo), dimension(atmos_auxvars2d_num), public atmos_auxvars2d_vinfo
integer, parameter, public atmos_diagvars_dens_id = 1
integer, parameter, public atmos_diagvars_u_id = 2
integer, parameter, public atmos_diagvars_v_id = 3
integer, parameter, public atmos_diagvars_w_id = 4
integer, parameter, public atmos_diagvars_t_id = 5
integer, parameter, public atmos_diagvars_umet_id = 6
integer, parameter, public atmos_diagvars_vmet_id = 7
integer, parameter, public atmos_diagvars_qdry_id = 8
integer, parameter, public atmos_diagvars_rh_id = 9
integer, parameter, public atmos_diagvars_engt_id = 10
integer, parameter, public atmos_diagvars_engp_id = 11
integer, parameter, public atmos_diagvars_engk_id = 12
integer, parameter, public atmos_diagvars_engi_id = 13
integer, parameter, public atmos_diagvars3d_num = 13
type(variableinfo), dimension(atmos_diagvars3d_num), public atmos_diagvars3d_vinfo
integer, parameter, public atmos_diagvars_rain_id = 1
integer, parameter, public atmos_diagvars_snow_id = 2
integer, parameter, public atmos_diagvars2d_num = 2
type(variableinfo), dimension(atmos_diagvars2d_num), public atmos_diagvars2d_vinfo

Detailed Description

module Atmosphere / Physics / Variables

Description
Container for atmospheric variables
Author
Yuta Kawai, Team SCALE
NAMELIST
  • PARAM_ATMOS_VARS
    nametypedefault valuecomment
    CHECK_RANGElogical.false.Flag whether the range of values is checked
    CHECK_TOTALlogical.false.

  • PARAM_ATMOS_VARS_RESTART
    nametypedefault valuecomment
    IN_BASENAMEcharacter(len=H_LONG)''Basename of the input file
    IN_POSTFIX_TIMELABELlogical.false.Add timelabel to the basename of input file?
    OUT_BASENAMEcharacter(len=H_LONG)''Basename of the output file
    OUT_POSTFIX_TIMELABELlogical.true.Add timelabel to the basename of output file?
    OUT_TITLEcharacter(len=H_MID)''Title of the output file
    OUT_DTYPEcharacter(len=H_SHORT)'DEFAULT'REAL4 or REAL8

History Output
namedescriptionunitvariable
PRECsurface precipitaion fluxkg/m2/sPREC
PREC_ENGIinternal energy of precipitationJ/m2PREC_ENGI
RAINsurface rain fluxkg/m2/sRAIN
SNOWsurface snow fluxkg/m2/sSNOW
DENSdensitykg/m3DENS
Uvelocity um/sU
Vvelocity vm/sV
Wvelocity wm/sW
TtemperatureKT
Umeteastward velocitym/sUmet
Vmetnorthward velocitym/sVmet

Function/Subroutine Documentation

◆ atmosvars_getlocalmeshprgvar()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshprgvar ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) prgvars_list,
class(modelvarmanager), intent(inout) auxvars_list,
integer, intent(in) varid,
class(localmeshfieldbase), intent(out), pointer var,
class(localmeshfieldbase), intent(out), optional, pointer dens_hyd,
class(localmeshfieldbase), intent(out), optional, pointer pres_hyd,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 898 of file mod_atmos_vars.F90.

901
902 implicit none
903
904 integer, intent(in) :: domID
905 class(MeshBase), intent(in) :: mesh
906 class(ModelVarManager), intent(inout) :: prgvars_list
907 class(ModelVarManager), intent(inout) :: auxvars_list
908 integer, intent(in) :: varid
909 class(LocalMeshFieldBase), pointer, intent(out) :: var
910 class(LocalMeshFieldBase), pointer, intent(out), optional :: DENS_hyd, PRES_hyd
911 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
912
913 class(MeshFieldBase), pointer :: field
914 class(LocalMeshBase), pointer :: lcmesh
915 !-------------------------------------------------------
916
917 !--
918 call prgvars_list%Get(varid, field)
919 call field%GetLocalMeshField(domid, var)
920
921 if (present(dens_hyd)) then
922 call auxvars_list%Get(auxvar_denshydro_id, field)
923 call field%GetLocalMeshField(domid, dens_hyd)
924 end if
925 if (present(pres_hyd)) then
926 call auxvars_list%Get(auxvar_preshydro_id, field)
927 call field%GetLocalMeshField(domid, pres_hyd)
928 end if
929
930 if (present(lcmesh3d)) then
931 call mesh%GetLocalMesh( domid, lcmesh )
932 nullify( lcmesh3d )
933
934 select type(lcmesh)
935 type is (localmesh3d)
936 if (present(lcmesh3d)) lcmesh3d => lcmesh
937 end select
938 end if
939
940 return

Referenced by mod_atmos_phy_tb::atmosphytb_setup(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshprgvars()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshprgvars ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) prgvars_list,
class(modelvarmanager), intent(inout) auxvars_list,
class(localmeshfieldbase), intent(out), pointer ddens,
class(localmeshfieldbase), intent(out), pointer momx,
class(localmeshfieldbase), intent(out), pointer momy,
class(localmeshfieldbase), intent(out), pointer momz,
class(localmeshfieldbase), intent(out), pointer therm,
class(localmeshfieldbase), intent(out), pointer dens_hyd,
class(localmeshfieldbase), intent(out), pointer pres_hyd,
class(localmeshfieldbase), intent(out), pointer rtot,
class(localmeshfieldbase), intent(out), pointer cvtot,
class(localmeshfieldbase), intent(out), pointer cptot,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 944 of file mod_atmos_vars.F90.

948 implicit none
949 integer, intent(in) :: domID
950 class(MeshBase), intent(in) :: mesh
951 class(ModelVarManager), intent(inout) :: prgvars_list
952 class(ModelVarManager), intent(inout) :: auxvars_list
953 class(LocalMeshFieldBase), pointer, intent(out) :: DDENS, MOMX, MOMY, MOMZ, THERM
954 class(LocalMeshFieldBase), pointer, intent(out) :: DENS_hyd, PRES_hyd
955 class(LocalMeshFieldBase), pointer, intent(out) :: Rtot, CVtot, CPtot
956 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
957
958 class(MeshFieldBase), pointer :: field
959 class(LocalMeshBase), pointer :: lcmesh
960 !-------------------------------------------------------
961
962 !--
963 call prgvars_list%Get(prgvar_ddens_id, field)
964 call field%GetLocalMeshField(domid, ddens)
965
966 call prgvars_list%Get(prgvar_momx_id, field)
967 call field%GetLocalMeshField(domid, momx)
968
969 call prgvars_list%Get(prgvar_momy_id, field)
970 call field%GetLocalMeshField(domid, momy)
971
972 call prgvars_list%Get(prgvar_momz_id, field)
973 call field%GetLocalMeshField(domid, momz)
974
975 call prgvars_list%Get(prgvar_therm_id, field)
976 call field%GetLocalMeshField(domid, therm)
977
978 !--
979 call auxvars_list%Get(auxvar_denshydro_id, field)
980 call field%GetLocalMeshField(domid, dens_hyd)
981
982 call auxvars_list%Get(auxvar_preshydro_id, field)
983 call field%GetLocalMeshField(domid, pres_hyd)
984
985 call auxvars_list%Get(auxvar_rtot_id, field)
986 call field%GetLocalMeshField(domid, rtot)
987
988 call auxvars_list%Get(auxvar_cvtot_id, field)
989 call field%GetLocalMeshField(domid, cvtot)
990
991 call auxvars_list%Get(auxvar_cptot_id, field)
992 call field%GetLocalMeshField(domid, cptot)
993
994 !---
995
996 if ( present(lcmesh3d) ) then
997 call mesh%GetLocalMesh( domid, lcmesh )
998 nullify( lcmesh3d )
999
1000 select type(lcmesh)
1001 type is (localmesh3d)
1002 if (present(lcmesh3d)) lcmesh3d => lcmesh
1003 end select
1004 end if
1005
1006 return

Referenced by mod_atmos_phy_tb::atmosphytb_setup(), mod_dg_prep::dg_prep(), mod_experiment::exp_setinitcond_lc::exp_setinitcond_lc(), mod_mkinit::mkinit(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshsfcvar()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshsfcvar ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) auxvars2d_list,
class(localmeshfieldbase), intent(out), pointer prec,
class(localmeshfieldbase), intent(out), pointer prec_engi,
class(localmesh2d), intent(out), optional, pointer lcmesh2d )

Definition at line 1010 of file mod_atmos_vars.F90.

1012
1013 implicit none
1014 integer, intent(in) :: domID
1015 class(MeshBase), intent(in) :: mesh
1016 class(ModelVarManager), intent(inout) :: auxvars2D_list
1017 class(LocalMeshFieldBase), pointer, intent(out) :: PREC, PREC_ENGI
1018 class(LocalMesh2D), pointer, intent(out), optional :: lcmesh2D
1019
1020 class(MeshFieldBase), pointer :: field
1021 class(LocalMeshBase), pointer :: lcmesh
1022 !-------------------------------------------------------
1023
1024 !--
1025 call auxvars2d_list%Get(atmos_auxvars2d_prec_id, field)
1026 call field%GetLocalMeshField(domid, prec)
1027
1028 call auxvars2d_list%Get(atmos_auxvars2d_prec_engi_id, field)
1029 call field%GetLocalMeshField(domid, prec_engi)
1030
1031 if (present(lcmesh2d)) then
1032 call mesh%GetLocalMesh( domid, lcmesh )
1033 nullify( lcmesh2d )
1034
1035 select type(lcmesh)
1036 type is (localmesh2d)
1037 if (present(lcmesh2d)) lcmesh2d => lcmesh
1038 end select
1039 end if
1040
1041 return

References atmos_auxvars2d_prec_engi_id, and atmos_auxvars2d_prec_id.

Referenced by mod_atmos_component::atmos_setup(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshqtrcvar()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshqtrcvar ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) trcvars_list,
integer, intent(in) varid,
class(localmeshfieldbase), intent(out), pointer var,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 1045 of file mod_atmos_vars.F90.

1048
1049 implicit none
1050 integer, intent(in) :: domID
1051 class(MeshBase), intent(in) :: mesh
1052 class(ModelVarManager), intent(inout) :: trcvars_list
1053 integer, intent(in) :: varid
1054 class(LocalMeshFieldBase), pointer, intent(out) :: var
1055 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1056
1057 class(MeshFieldBase), pointer :: field
1058 class(LocalMeshBase), pointer :: lcmesh
1059 !-------------------------------------------------------
1060
1061 !--
1062 call trcvars_list%Get(varid, field)
1063 call field%GetLocalMeshField(domid, var)
1064
1065 if (present(lcmesh3d)) then
1066 call mesh%GetLocalMesh( domid, lcmesh )
1067 nullify( lcmesh3d )
1068
1069 select type(lcmesh)
1070 type is (localmesh3d)
1071 if (present(lcmesh3d)) lcmesh3d => lcmesh
1072 end select
1073 end if
1074
1075 return

Referenced by mod_atmos_phy_tb::atmosphytb_setup(), mod_experiment::exp_setinitcond_lc::exp_setinitcond_lc(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshqtrc_qv()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshqtrc_qv ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) trcvars_list,
class(modelvarmanager), intent(inout) forcing_list,
class(localmeshfieldbase), intent(out), pointer var,
class(localmeshfieldbase), intent(out), pointer var_tp,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 1079 of file mod_atmos_vars.F90.

1081
1082 use scale_atmos_hydrometeor, only: &
1083 atmos_hydrometeor_dry, &
1084 i_qv
1085 implicit none
1086 integer, intent(in) :: domID
1087 class(MeshBase), intent(in) :: mesh
1088 class(ModelVarManager), intent(inout) :: trcvars_list
1089 class(ModelVarManager), intent(inout) :: forcing_list
1090 class(LocalMeshFieldBase), pointer, intent(out) :: var
1091 class(LocalMeshFieldBase), pointer, intent(out) :: var_tp
1092 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1093
1094 class(MeshFieldBase), pointer :: field
1095 class(LocalMeshBase), pointer :: lcmesh
1096
1097 integer :: iq, tend_iq
1098 !-------------------------------------------------------
1099
1100 !--
1101
1102 if ( atmos_hydrometeor_dry ) then
1103 iq = 0; tend_iq = phytend_num1+1
1104 else
1105 iq = i_qv; tend_iq = phytend_num1 + i_qv
1106 end if
1107
1108 call trcvars_list%Get(iq, field)
1109 call field%GetLocalMeshField(domid, var)
1110
1111 call forcing_list%Get(tend_iq, field)
1112 call field%GetLocalMeshField(domid, var_tp)
1113
1114 if (present(lcmesh3d)) then
1115 call mesh%GetLocalMesh( domid, lcmesh )
1116 nullify( lcmesh3d )
1117
1118 select type(lcmesh)
1119 type is (localmesh3d)
1120 if (present(lcmesh3d)) lcmesh3d => lcmesh
1121 end select
1122 end if
1123
1124 return

Referenced by mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshqtrcvarlist()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshqtrcvarlist ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) trcvars_list,
integer, intent(in) varid_s,
type(localmeshfieldbaselist), dimension(:), intent(out) var_list,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 1128 of file mod_atmos_vars.F90.

1131
1132 implicit none
1133 integer, intent(in) :: domID
1134 class(MeshBase), intent(in) :: mesh
1135 class(ModelVarManager), intent(inout) :: trcvars_list
1136 integer, intent(in) :: varid_s
1137 type(LocalMeshFieldBaseList), intent(out) :: var_list(:)
1138 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1139
1140 class(MeshFieldBase), pointer :: field
1141 class(LocalMeshBase), pointer :: lcmesh
1142
1143 integer :: iq
1144 !-------------------------------------------------------
1145
1146 !--
1147 do iq = varid_s, varid_s + size(var_list) - 1
1148 call trcvars_list%Get(iq, field)
1149 call field%GetLocalMeshField(domid, var_list(iq-varid_s+1)%ptr)
1150 end do
1151 if (present(lcmesh3d)) then
1152 call mesh%GetLocalMesh( domid, lcmesh )
1153 nullify( lcmesh3d )
1154
1155 select type(lcmesh)
1156 type is (localmesh3d)
1157 if (present(lcmesh3d)) lcmesh3d => lcmesh
1158 end select
1159 end if
1160
1161 return

Referenced by atmosvars_getlocalmeshqtrcphytend(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshphyauxvars()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshphyauxvars ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) phyauxvars_list,
class(localmeshfieldbase), intent(out), pointer pres,
class(localmeshfieldbase), intent(out), pointer pt,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 1165 of file mod_atmos_vars.F90.

1168
1169 implicit none
1170 integer, intent(in) :: domID
1171 class(MeshBase), intent(in) :: mesh
1172 class(ModelVarManager), intent(inout) :: phyauxvars_list
1173 class(LocalMeshFieldBase), pointer, intent(out) :: PRES, PT
1174 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1175
1176 class(MeshFieldBase), pointer :: field
1177 class(LocalMeshBase), pointer :: lcmesh
1178 !-------------------------------------------------------
1179
1180 !--
1181 call phyauxvars_list%Get(auxvar_pres_id, field)
1182 call field%GetLocalMeshField(domid, pres)
1183
1184 call phyauxvars_list%Get(auxvar_pt_id, field)
1185 call field%GetLocalMeshField(domid, pt)
1186
1187 !---
1188
1189 if (present(lcmesh3d)) then
1190 call mesh%GetLocalMesh( domid, lcmesh )
1191 nullify( lcmesh3d )
1192
1193 select type(lcmesh)
1194 type is (localmesh3d)
1195 if (present(lcmesh3d)) lcmesh3d => lcmesh
1196 end select
1197 end if
1198
1199 return

Referenced by mod_atmos_phy_tb::atmosphytb_setup(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshphytends()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshphytends ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) phytends_list,
class(localmeshfieldbase), intent(out), pointer dens_tp,
class(localmeshfieldbase), intent(out), pointer momx_tp,
class(localmeshfieldbase), intent(out), pointer momy_tp,
class(localmeshfieldbase), intent(out), pointer momz_tp,
class(localmeshfieldbase), intent(out), pointer rhot_tp,
class(localmeshfieldbase), intent(out), pointer rhoh_p,
type(localmeshfieldbaselist), dimension(qa), intent(inout), optional rhoq_tp,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 1203 of file mod_atmos_vars.F90.

1207
1208 implicit none
1209 integer, intent(in) :: domID
1210 class(MeshBase), intent(in) :: mesh
1211 class(ModelVarManager), intent(inout) :: phytends_list
1212 class(LocalMeshFieldBase), pointer, intent(out) :: DENS_tp, MOMX_tp, MOMY_tp, MOMZ_tp, RHOT_tp
1213 class(LocalMeshFieldBase), pointer, intent(out) :: RHOH_p
1214 type(LocalMeshFieldBaseList), intent(inout), optional :: RHOQ_tp(QA)
1215 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1216
1217 class(MeshFieldBase), pointer :: field
1218 class(LocalMeshBase), pointer :: lcmesh
1219
1220 integer :: iq
1221 !-------------------------------------------------------
1222
1223 !--
1224 call phytends_list%Get(phytend_dens_id, field)
1225 call field%GetLocalMeshField(domid, dens_tp)
1226
1227 call phytends_list%Get(phytend_momx_id, field)
1228 call field%GetLocalMeshField(domid, momx_tp)
1229
1230 call phytends_list%Get(phytend_momy_id, field)
1231 call field%GetLocalMeshField(domid, momy_tp)
1232
1233 call phytends_list%Get(phytend_momz_id, field)
1234 call field%GetLocalMeshField(domid, momz_tp)
1235
1236 call phytends_list%Get(phytend_rhot_id, field)
1237 call field%GetLocalMeshField(domid, rhot_tp)
1238
1239 call phytends_list%Get(phytend_rhoh_id, field)
1240 call field%GetLocalMeshField(domid, rhoh_p)
1241
1242 if ( present(rhoq_tp) ) then
1243 do iq = 1, qa
1244 call phytends_list%Get(phytend_num1+iq, field)
1245 call field%GetLocalMeshField(domid, rhoq_tp(iq)%ptr)
1246 end do
1247 end if
1248
1249 !---
1250 if ( present(lcmesh3d) ) then
1251 call mesh%GetLocalMesh( domid, lcmesh )
1252 nullify( lcmesh3d )
1253
1254 select type(lcmesh)
1255 type is (localmesh3d)
1256 if ( present(lcmesh3d) ) lcmesh3d => lcmesh
1257 end select
1258 end if
1259
1260 return

Referenced by mod_atmos_component::atmos_setup(), mod_atmos_phy_tb::atmosphytb_setup(), and mod_atmos_vars::atmosvars::regist_physvar_manager().

◆ atmosvars_getlocalmeshqtrcphytend()

subroutine, public mod_atmos_vars::atmosvars_getlocalmeshqtrcphytend ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) phytends_list,
integer, intent(in) qtrcid,
class(localmeshfieldbase), intent(out), pointer rhoq_tp )

Definition at line 1264 of file mod_atmos_vars.F90.

1267
1268 implicit none
1269 integer, intent(in) :: domID
1270 class(MeshBase), intent(in) :: mesh
1271 class(ModelVarManager), intent(inout) :: phytends_list
1272 integer, intent(in) :: qtrcid
1273 class(LocalMeshFieldBase), pointer, intent(out) :: RHOQ_tp
1274
1275 class(MeshFieldBase), pointer :: field
1276 class(LocalMeshBase), pointer :: lcmesh
1277 !-------------------------------------------------------
1278
1279 call phytends_list%Get(phytend_num1 + qtrcid, field)
1280 call field%GetLocalMeshField(domid, rhoq_tp)
1281
1282 return

References mod_atmos_phy_mp_vars::atmosphympvars_getlocalmeshfields_sfcflx(), and atmosvars_getlocalmeshqtrcvarlist().

Referenced by mod_atmos_vars::atmosvars::regist_physvar_manager().

Variable Documentation

◆ progvars_check_min

real(rp), dimension(prgvar_num), parameter mod_atmos_vars::progvars_check_min = (/ -1.0_RP, -100.0_RP, -200.0_RP, -200.0_RP, -200.0_RP /)

Definition at line 144 of file mod_atmos_vars.F90.

144 real(RP), parameter :: PROGVARS_check_min(PRGVAR_NUM) = (/ -1.0_rp, -100.0_rp, -200.0_rp, -200.0_rp, -200.0_rp /)

◆ atmos_auxvars2d_prec_id

integer, parameter, public mod_atmos_vars::atmos_auxvars2d_prec_id = 1

Definition at line 150 of file mod_atmos_vars.F90.

150 integer, public, parameter :: ATMOS_AUXVARS2D_PREC_ID = 1

Referenced by atmosvars_getlocalmeshsfcvar().

◆ atmos_auxvars2d_prec_engi_id

integer, parameter, public mod_atmos_vars::atmos_auxvars2d_prec_engi_id = 2

Definition at line 151 of file mod_atmos_vars.F90.

151 integer, public, parameter :: ATMOS_AUXVARS2D_PREC_ENGI_ID = 2

Referenced by atmosvars_getlocalmeshsfcvar().

◆ atmos_auxvars2d_num

integer, parameter, public mod_atmos_vars::atmos_auxvars2d_num = 2

Definition at line 152 of file mod_atmos_vars.F90.

152 integer, public, parameter :: ATMOS_AUXVARS2D_NUM = 2

◆ atmos_auxvars2d_vinfo

type(variableinfo), dimension(atmos_auxvars2d_num), public mod_atmos_vars::atmos_auxvars2d_vinfo

Definition at line 154 of file mod_atmos_vars.F90.

154 type(VariableInfo), public :: ATMOS_AUXVARS2D_VINFO(ATMOS_AUXVARS2D_NUM)

◆ atmos_diagvars_dens_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_dens_id = 1

Definition at line 163 of file mod_atmos_vars.F90.

163 integer, public, parameter :: ATMOS_DIAGVARS_DENS_ID = 1

◆ atmos_diagvars_u_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_u_id = 2

Definition at line 164 of file mod_atmos_vars.F90.

164 integer, public, parameter :: ATMOS_DIAGVARS_U_ID = 2

◆ atmos_diagvars_v_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_v_id = 3

Definition at line 165 of file mod_atmos_vars.F90.

165 integer, public, parameter :: ATMOS_DIAGVARS_V_ID = 3

◆ atmos_diagvars_w_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_w_id = 4

Definition at line 166 of file mod_atmos_vars.F90.

166 integer, public, parameter :: ATMOS_DIAGVARS_W_ID = 4

◆ atmos_diagvars_t_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_t_id = 5

Definition at line 167 of file mod_atmos_vars.F90.

167 integer, public, parameter :: ATMOS_DIAGVARS_T_ID = 5

◆ atmos_diagvars_umet_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_umet_id = 6

Definition at line 168 of file mod_atmos_vars.F90.

168 integer, public, parameter :: ATMOS_DIAGVARS_Umet_ID = 6

◆ atmos_diagvars_vmet_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_vmet_id = 7

Definition at line 169 of file mod_atmos_vars.F90.

169 integer, public, parameter :: ATMOS_DIAGVARS_Vmet_ID = 7

◆ atmos_diagvars_qdry_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_qdry_id = 8

Definition at line 170 of file mod_atmos_vars.F90.

170 integer, public, parameter :: ATMOS_DIAGVARS_QDRY_ID = 8

◆ atmos_diagvars_rh_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_rh_id = 9

Definition at line 171 of file mod_atmos_vars.F90.

171 integer, public, parameter :: ATMOS_DIAGVARS_RH_ID = 9

◆ atmos_diagvars_engt_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_engt_id = 10

Definition at line 172 of file mod_atmos_vars.F90.

172 integer, public, parameter :: ATMOS_DIAGVARS_ENGT_ID = 10

◆ atmos_diagvars_engp_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_engp_id = 11

Definition at line 173 of file mod_atmos_vars.F90.

173 integer, public, parameter :: ATMOS_DIAGVARS_ENGP_ID = 11

◆ atmos_diagvars_engk_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_engk_id = 12

Definition at line 174 of file mod_atmos_vars.F90.

174 integer, public, parameter :: ATMOS_DIAGVARS_ENGK_ID = 12

◆ atmos_diagvars_engi_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_engi_id = 13

Definition at line 175 of file mod_atmos_vars.F90.

175 integer, public, parameter :: ATMOS_DIAGVARS_ENGI_ID = 13

◆ atmos_diagvars3d_num

integer, parameter, public mod_atmos_vars::atmos_diagvars3d_num = 13

Definition at line 176 of file mod_atmos_vars.F90.

176 integer, public, parameter :: ATMOS_DIAGVARS3D_NUM = 13

◆ atmos_diagvars3d_vinfo

type(variableinfo), dimension(atmos_diagvars3d_num), public mod_atmos_vars::atmos_diagvars3d_vinfo

Definition at line 178 of file mod_atmos_vars.F90.

178 type(VariableInfo), public :: ATMOS_DIAGVARS3D_VINFO(ATMOS_DIAGVARS3D_NUM)

◆ atmos_diagvars_rain_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_rain_id = 1

Definition at line 195 of file mod_atmos_vars.F90.

195 integer, public, parameter :: ATMOS_DIAGVARS_RAIN_ID = 1

◆ atmos_diagvars_snow_id

integer, parameter, public mod_atmos_vars::atmos_diagvars_snow_id = 2

Definition at line 196 of file mod_atmos_vars.F90.

196 integer, public, parameter :: ATMOS_DIAGVARS_SNOW_ID = 2

◆ atmos_diagvars2d_num

integer, parameter, public mod_atmos_vars::atmos_diagvars2d_num = 2

Definition at line 197 of file mod_atmos_vars.F90.

197 integer, public, parameter :: ATMOS_DIAGVARS2D_NUM = 2

◆ atmos_diagvars2d_vinfo

type(variableinfo), dimension(atmos_diagvars2d_num), public mod_atmos_vars::atmos_diagvars2d_vinfo

Definition at line 199 of file mod_atmos_vars.F90.

199 type(VariableInfo), public :: ATMOS_DIAGVARS2D_VINFO(ATMOS_DIAGVARS2D_NUM)