FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
mod_atmos_vars Module Reference

module ATMOSPHERE / Variables More...

Data Types

type  atmosvars
 

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 / Variables

Description
Container for atmospheric variables
Author
Yuta Kawai, Team SCALE
NAMELIST
  • PARAM_ATMOS_VARS
    nametypedefault valuecomment
    CHECK_RANGE logical .false.
    CHECK_TOTAL logical .false.

  • PARAM_ATMOS_VARS_RESTART
    nametypedefault valuecomment
    IN_BASENAME character(len=H_LONG) '' Basename of the input file
    IN_POSTFIX_TIMELABEL logical .false. Add timelabel to the basename of input file?
    OUT_BASENAME character(len=H_LONG) '' Basename of the output file
    OUT_POSTFIX_TIMELABEL logical .true. Add timelabel to the basename of output file?
    OUT_TITLE character(len=H_MID) '' Title of the output file
    OUT_DTYPE character(len=H_SHORT) 'DEFAULT' REAL4 or REAL8

History Output
namedescriptionunitvariable
PREC surface precipitaion flux kg/m2/s PREC
PREC_ENGI internal energy of precipitation J/m2 PREC_ENGI
RAIN surface rain flux kg/m2/s RAIN
SNOW surface snow flux kg/m2/s SNOW
DENS density kg/m3 DENS
U velocity u m/s U
V velocity v m/s V
W velocity w m/s W
T temperature K T
Umet eastward velocity m/s Umet
Vmet northward velocity m/s Vmet

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 873 of file mod_atmos_vars.F90.

876
877 implicit none
878
879 integer, intent(in) :: domID
880 class(MeshBase), intent(in) :: mesh
881 class(ModelVarManager), intent(inout) :: prgvars_list
882 class(ModelVarManager), intent(inout) :: auxvars_list
883 integer, intent(in) :: varid
884 class(LocalMeshFieldBase), pointer, intent(out) :: var
885 class(LocalMeshFieldBase), pointer, intent(out), optional :: DENS_hyd, PRES_hyd
886 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
887
888 class(MeshFieldBase), pointer :: field
889 class(LocalMeshBase), pointer :: lcmesh
890 !-------------------------------------------------------
891
892 !--
893 call prgvars_list%Get(varid, field)
894 call field%GetLocalMeshField(domid, var)
895
896 if (present(dens_hyd)) then
897 call auxvars_list%Get(auxvar_denshydro_id, field)
898 call field%GetLocalMeshField(domid, dens_hyd)
899 end if
900 if (present(pres_hyd)) then
901 call auxvars_list%Get(auxvar_preshydro_id, field)
902 call field%GetLocalMeshField(domid, pres_hyd)
903 end if
904
905 if (present(lcmesh3d)) then
906 call mesh%GetLocalMesh( domid, lcmesh )
907 nullify( lcmesh3d )
908
909 select type(lcmesh)
910 type is (localmesh3d)
911 if (present(lcmesh3d)) lcmesh3d => lcmesh
912 end select
913 end if
914
915 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 919 of file mod_atmos_vars.F90.

923 implicit none
924 integer, intent(in) :: domID
925 class(MeshBase), intent(in) :: mesh
926 class(ModelVarManager), intent(inout) :: prgvars_list
927 class(ModelVarManager), intent(inout) :: auxvars_list
928 class(LocalMeshFieldBase), pointer, intent(out) :: DDENS, MOMX, MOMY, MOMZ, THERM
929 class(LocalMeshFieldBase), pointer, intent(out) :: DENS_hyd, PRES_hyd
930 class(LocalMeshFieldBase), pointer, intent(out) :: Rtot, CVtot, CPtot
931 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
932
933 class(MeshFieldBase), pointer :: field
934 class(LocalMeshBase), pointer :: lcmesh
935 !-------------------------------------------------------
936
937 !--
938 call prgvars_list%Get(prgvar_ddens_id, field)
939 call field%GetLocalMeshField(domid, ddens)
940
941 call prgvars_list%Get(prgvar_momx_id, field)
942 call field%GetLocalMeshField(domid, momx)
943
944 call prgvars_list%Get(prgvar_momy_id, field)
945 call field%GetLocalMeshField(domid, momy)
946
947 call prgvars_list%Get(prgvar_momz_id, field)
948 call field%GetLocalMeshField(domid, momz)
949
950 call prgvars_list%Get(prgvar_therm_id, field)
951 call field%GetLocalMeshField(domid, therm)
952
953 !--
954 call auxvars_list%Get(auxvar_denshydro_id, field)
955 call field%GetLocalMeshField(domid, dens_hyd)
956
957 call auxvars_list%Get(auxvar_preshydro_id, field)
958 call field%GetLocalMeshField(domid, pres_hyd)
959
960 call auxvars_list%Get(auxvar_rtot_id, field)
961 call field%GetLocalMeshField(domid, rtot)
962
963 call auxvars_list%Get(auxvar_cvtot_id, field)
964 call field%GetLocalMeshField(domid, cvtot)
965
966 call auxvars_list%Get(auxvar_cptot_id, field)
967 call field%GetLocalMeshField(domid, cptot)
968
969 !---
970
971 if ( present(lcmesh3d) ) then
972 call mesh%GetLocalMesh( domid, lcmesh )
973 nullify( lcmesh3d )
974
975 select type(lcmesh)
976 type is (localmesh3d)
977 if (present(lcmesh3d)) lcmesh3d => lcmesh
978 end select
979 end if
980
981 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 985 of file mod_atmos_vars.F90.

987
988 implicit none
989 integer, intent(in) :: domID
990 class(MeshBase), intent(in) :: mesh
991 class(ModelVarManager), intent(inout) :: auxvars2D_list
992 class(LocalMeshFieldBase), pointer, intent(out) :: PREC, PREC_ENGI
993 class(LocalMesh2D), pointer, intent(out), optional :: lcmesh2D
994
995 class(MeshFieldBase), pointer :: field
996 class(LocalMeshBase), pointer :: lcmesh
997 !-------------------------------------------------------
998
999 !--
1000 call auxvars2d_list%Get(atmos_auxvars2d_prec_id, field)
1001 call field%GetLocalMeshField(domid, prec)
1002
1003 call auxvars2d_list%Get(atmos_auxvars2d_prec_engi_id, field)
1004 call field%GetLocalMeshField(domid, prec_engi)
1005
1006 if (present(lcmesh2d)) then
1007 call mesh%GetLocalMesh( domid, lcmesh )
1008 nullify( lcmesh2d )
1009
1010 select type(lcmesh)
1011 type is (localmesh2d)
1012 if (present(lcmesh2d)) lcmesh2d => lcmesh
1013 end select
1014 end if
1015
1016 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 1020 of file mod_atmos_vars.F90.

1023
1024 implicit none
1025 integer, intent(in) :: domID
1026 class(MeshBase), intent(in) :: mesh
1027 class(ModelVarManager), intent(inout) :: trcvars_list
1028 integer, intent(in) :: varid
1029 class(LocalMeshFieldBase), pointer, intent(out) :: var
1030 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1031
1032 class(MeshFieldBase), pointer :: field
1033 class(LocalMeshBase), pointer :: lcmesh
1034 !-------------------------------------------------------
1035
1036 !--
1037 call trcvars_list%Get(varid, field)
1038 call field%GetLocalMeshField(domid, var)
1039
1040 if (present(lcmesh3d)) then
1041 call mesh%GetLocalMesh( domid, lcmesh )
1042 nullify( lcmesh3d )
1043
1044 select type(lcmesh)
1045 type is (localmesh3d)
1046 if (present(lcmesh3d)) lcmesh3d => lcmesh
1047 end select
1048 end if
1049
1050 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 1054 of file mod_atmos_vars.F90.

1056
1057 use scale_atmos_hydrometeor, only: &
1058 atmos_hydrometeor_dry, &
1059 i_qv
1060 implicit none
1061 integer, intent(in) :: domID
1062 class(MeshBase), intent(in) :: mesh
1063 class(ModelVarManager), intent(inout) :: trcvars_list
1064 class(ModelVarManager), intent(inout) :: forcing_list
1065 class(LocalMeshFieldBase), pointer, intent(out) :: var
1066 class(LocalMeshFieldBase), pointer, intent(out) :: var_tp
1067 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1068
1069 class(MeshFieldBase), pointer :: field
1070 class(LocalMeshBase), pointer :: lcmesh
1071
1072 integer :: iq, tend_iq
1073 !-------------------------------------------------------
1074
1075 !--
1076 if ( atmos_hydrometeor_dry ) then
1077 iq =0; tend_iq = phytend_num1+1
1078 else
1079 iq = i_qv; tend_iq = phytend_num1 + i_qv
1080 end if
1081
1082 call trcvars_list%Get(iq, field)
1083 call field%GetLocalMeshField(domid, var)
1084
1085 call forcing_list%Get(tend_iq, field)
1086 call field%GetLocalMeshField(domid, var_tp)
1087
1088 if (present(lcmesh3d)) then
1089 call mesh%GetLocalMesh( domid, lcmesh )
1090 nullify( lcmesh3d )
1091
1092 select type(lcmesh)
1093 type is (localmesh3d)
1094 if (present(lcmesh3d)) lcmesh3d => lcmesh
1095 end select
1096 end if
1097
1098 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 1102 of file mod_atmos_vars.F90.

1105
1106 implicit none
1107 integer, intent(in) :: domID
1108 class(MeshBase), intent(in) :: mesh
1109 class(ModelVarManager), intent(inout) :: trcvars_list
1110 integer, intent(in) :: varid_s
1111 type(LocalMeshFieldBaseList), intent(out) :: var_list(:)
1112 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1113
1114 class(MeshFieldBase), pointer :: field
1115 class(LocalMeshBase), pointer :: lcmesh
1116
1117 integer :: iq
1118 !-------------------------------------------------------
1119
1120 !--
1121 do iq = varid_s, varid_s + size(var_list) - 1
1122 call trcvars_list%Get(iq, field)
1123 call field%GetLocalMeshField(domid, var_list(iq-varid_s+1)%ptr)
1124 end do
1125 if (present(lcmesh3d)) then
1126 call mesh%GetLocalMesh( domid, lcmesh )
1127 nullify( lcmesh3d )
1128
1129 select type(lcmesh)
1130 type is (localmesh3d)
1131 if (present(lcmesh3d)) lcmesh3d => lcmesh
1132 end select
1133 end if
1134
1135 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 1139 of file mod_atmos_vars.F90.

1142
1143 implicit none
1144 integer, intent(in) :: domID
1145 class(MeshBase), intent(in) :: mesh
1146 class(ModelVarManager), intent(inout) :: phyauxvars_list
1147 class(LocalMeshFieldBase), pointer, intent(out) :: PRES, PT
1148 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1149
1150 class(MeshFieldBase), pointer :: field
1151 class(LocalMeshBase), pointer :: lcmesh
1152 !-------------------------------------------------------
1153
1154 !--
1155 call phyauxvars_list%Get(auxvar_pres_id, field)
1156 call field%GetLocalMeshField(domid, pres)
1157
1158 call phyauxvars_list%Get(auxvar_pt_id, field)
1159 call field%GetLocalMeshField(domid, pt)
1160
1161 !---
1162
1163 if (present(lcmesh3d)) then
1164 call mesh%GetLocalMesh( domid, lcmesh )
1165 nullify( lcmesh3d )
1166
1167 select type(lcmesh)
1168 type is (localmesh3d)
1169 if (present(lcmesh3d)) lcmesh3d => lcmesh
1170 end select
1171 end if
1172
1173 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 1177 of file mod_atmos_vars.F90.

1181
1182 implicit none
1183 integer, intent(in) :: domID
1184 class(MeshBase), intent(in) :: mesh
1185 class(ModelVarManager), intent(inout) :: phytends_list
1186 class(LocalMeshFieldBase), pointer, intent(out) :: DENS_tp, MOMX_tp, MOMY_tp, MOMZ_tp, RHOT_tp
1187 class(LocalMeshFieldBase), pointer, intent(out) :: RHOH_p
1188 type(LocalMeshFieldBaseList), intent(inout), optional :: RHOQ_tp(QA)
1189 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
1190
1191 class(MeshFieldBase), pointer :: field
1192 class(LocalMeshBase), pointer :: lcmesh
1193
1194 integer :: iq
1195 !-------------------------------------------------------
1196
1197 !--
1198 call phytends_list%Get(phytend_dens_id, field)
1199 call field%GetLocalMeshField(domid, dens_tp)
1200
1201 call phytends_list%Get(phytend_momx_id, field)
1202 call field%GetLocalMeshField(domid, momx_tp)
1203
1204 call phytends_list%Get(phytend_momy_id, field)
1205 call field%GetLocalMeshField(domid, momy_tp)
1206
1207 call phytends_list%Get(phytend_momz_id, field)
1208 call field%GetLocalMeshField(domid, momz_tp)
1209
1210 call phytends_list%Get(phytend_rhot_id, field)
1211 call field%GetLocalMeshField(domid, rhot_tp)
1212
1213 call phytends_list%Get(phytend_rhoh_id, field)
1214 call field%GetLocalMeshField(domid, rhoh_p)
1215
1216 if ( present(rhoq_tp) ) then
1217 do iq = 1, qa
1218 call phytends_list%Get(phytend_num1+iq, field)
1219 call field%GetLocalMeshField(domid, rhoq_tp(iq)%ptr)
1220 end do
1221 end if
1222
1223 !---
1224 if ( present(lcmesh3d) ) then
1225 call mesh%GetLocalMesh( domid, lcmesh )
1226 nullify( lcmesh3d )
1227
1228 select type(lcmesh)
1229 type is (localmesh3d)
1230 if ( present(lcmesh3d) ) lcmesh3d => lcmesh
1231 end select
1232 end if
1233
1234 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 1238 of file mod_atmos_vars.F90.

1241
1242 implicit none
1243 integer, intent(in) :: domID
1244 class(MeshBase), intent(in) :: mesh
1245 class(ModelVarManager), intent(inout) :: phytends_list
1246 integer, intent(in) :: qtrcid
1247 class(LocalMeshFieldBase), pointer, intent(out) :: RHOQ_tp
1248
1249 class(MeshFieldBase), pointer :: field
1250 class(LocalMeshBase), pointer :: lcmesh
1251 !-------------------------------------------------------
1252
1253 call phytends_list%Get(phytend_num1 + qtrcid, field)
1254 call field%GetLocalMeshField(domid, rhoq_tp)
1255
1256 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 142 of file mod_atmos_vars.F90.

142 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 148 of file mod_atmos_vars.F90.

148 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 149 of file mod_atmos_vars.F90.

149 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 150 of file mod_atmos_vars.F90.

150 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 152 of file mod_atmos_vars.F90.

152 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 161 of file mod_atmos_vars.F90.

161 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 162 of file mod_atmos_vars.F90.

162 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 163 of file mod_atmos_vars.F90.

163 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 164 of file mod_atmos_vars.F90.

164 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 165 of file mod_atmos_vars.F90.

165 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 166 of file mod_atmos_vars.F90.

166 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 167 of file mod_atmos_vars.F90.

167 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 168 of file mod_atmos_vars.F90.

168 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 169 of file mod_atmos_vars.F90.

169 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 170 of file mod_atmos_vars.F90.

170 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 171 of file mod_atmos_vars.F90.

171 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 172 of file mod_atmos_vars.F90.

172 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 173 of file mod_atmos_vars.F90.

173 integer, public, parameter :: ATMOS_DIAGVARS_ENGI_ID = 13

◆ atmos_diagvars3d_num

integer, parameter, public mod_atmos_vars::atmos_diagvars3d_num = 13

Definition at line 174 of file mod_atmos_vars.F90.

174 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 176 of file mod_atmos_vars.F90.

176 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 193 of file mod_atmos_vars.F90.

193 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 194 of file mod_atmos_vars.F90.

194 integer, public, parameter :: ATMOS_DIAGVARS_SNOW_ID = 2

◆ atmos_diagvars2d_num

integer, parameter, public mod_atmos_vars::atmos_diagvars2d_num = 2

Definition at line 195 of file mod_atmos_vars.F90.

195 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 197 of file mod_atmos_vars.F90.

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