11#include "scaleFElib.h"
21 use scale_tracer,
only: &
22 qa, tracer_name, tracer_desc, tracer_unit
23 use scale_atmos_hydrometeor,
only: &
78 integer :: prog_vars_commid
83 integer :: qtrc_vars_commid
88 integer :: aux_vars_commid
104 integer :: phytends_commid
105 integer :: phytend_num_tot
108 integer,
allocatable :: diagvars2d_histid(:)
109 integer,
allocatable :: diagvars3d_histid(:)
113 logical :: check_range
114 logical :: check_total
117 procedure :: init => atmosvars_init
118 procedure :: final => atmosvars_final
119 procedure :: calc_diagnostics => atmosvars_calculatediagnostics
120 procedure :: calc_diagvar => atmosvars_calcdiagvar
121 procedure :: calc_diagvar2d => atmosvars_calcdiagvar2d
122 procedure :: history => atmosvars_history
123 procedure :: check => atmosvars_check
124 procedure :: monitor => atmosvars_monitor
125 procedure :: read_restart_file => atmosvar_read_restart_file
126 procedure :: write_restart_file => atmosvar_write_restart_file
127 procedure :: regist_physvar_manager => atmosvars_regist_physvar_manager
145 real(rp),
parameter :: progvars_check_max(
prgvar_num) = (/ 1.0_rp, 100.0_rp, 200.0_rp, 200.0_rp, 200.0_rp /)
156 variableinfo(
atmos_auxvars2d_prec_id ,
'PREC',
'surface precipitaion flux' ,
'kg/m2/s', 2,
'XY',
'precipitation_flux' ), &
183 variableinfo(
atmos_diagvars_w_id ,
'W' ,
'velocity w' ,
'm/s' , 3,
'XYZ',
'upward_air_velocity' ), &
184 variableinfo(
atmos_diagvars_t_id ,
'T' ,
'temperature' ,
'K' , 3,
'XYZ',
'air_temperature' ), &
188 variableinfo(
atmos_diagvars_rh_id ,
'RH',
'relative humidity(liq)',
'%', 3,
'XYZ',
'relative_humidity' ), &
201 variableinfo(
atmos_diagvars_rain_id,
'RAIN',
'surface rain flux' ,
'kg/m2/s', 2,
'XY',
'rainfall_flux' ), &
202 variableinfo(
atmos_diagvars_snow_id,
'SNOW',
'surface snow flux' ,
'kg/m2/s', 2,
'XY',
'snowfall_flux' ) /
210 private :: vars_calc_diagnosevar_lc
211 private :: vars_calc_diagnosevar2d_lc
215 integer,
private,
parameter :: im_qdry = 1
216 integer,
private,
parameter :: im_qtot = 2
217 integer,
private,
parameter :: im_engt = 3
218 integer,
private,
parameter :: im_engp = 4
219 integer,
private,
parameter :: im_engk = 5
220 integer,
private,
parameter :: im_engi = 6
221 integer,
private,
parameter :: dvm_nmax = 6
222 integer,
private :: dv_monit_id(dvm_nmax)
231 subroutine atmosvars_init( this, atm_mesh )
232 use scale_atmos_hydrometeor,
only: &
233 atmos_hydrometeor_dry
235 monitor_reg => file_monitor_meshfield_reg
242 class(
atmosvars),
target,
intent(inout) :: this
243 class(
atmosmesh),
target,
intent(inout) :: atm_mesh
248 logical :: reg_file_hist
254 logical :: check_range = .false.
255 logical :: check_total = .false.
257 namelist / param_atmos_vars / &
261 character(len=H_LONG) :: in_basename =
''
262 logical :: in_postfix_timelabel = .false.
263 character(len=H_LONG) :: out_basename =
''
264 logical :: out_postfix_timelabel = .true.
265 character(len=H_MID) :: out_title =
''
266 character(len=H_SHORT) :: out_dtype =
'DEFAULT'
268 namelist / param_atmos_vars_restart / &
270 in_postfix_timelabel, &
272 out_postfix_timelabel, &
277 logical :: is_specified
285 log_info(
'AtmosVars_Init',*)
289 read(io_fid_conf,nml=param_atmos_vars,iostat=ierr)
291 log_info(
"ATMOS_vars_setup",*)
'Not found namelist. Default used.'
292 elseif( ierr > 0 )
then
293 log_error(
"ATMOS_vars_setup",*)
'Invalid names in namelist PARAM_ATMOS_VARS. Check!'
296 log_nml(param_atmos_vars)
299 this%mesh => atm_mesh
300 mesh3d => atm_mesh%ptr_mesh
301 call mesh3d%GetMesh2D( mesh2d )
306 call this%PROGVARS_manager%Init()
307 call this%QTRCVARS_manager%Init()
308 call this%AUXVARS_manager%Init()
309 call this%PHYTENDS_manager%Init()
312 allocate( this%QTRC_VARS(0:qa) )
315 this%PHYTEND_NUM_TOT = phytend_num1 + max(1,qa)
316 allocate( this%PHY_TEND(this%PHYTEND_NUM_TOT) )
319 this%PROG_VARS, this%QTRC_VARS, this%AUX_VARS, this%PHY_TEND, &
320 this%PROGVARS_manager, this%QTRCVARS_manager, this%AUXVARS_manager, this%PHYTENDS_manager, &
321 this%PHYTEND_NUM_TOT, mesh3d, &
326 call atm_mesh%Create_communicator( &
328 this%PROGVARS_manager, &
330 this%PROG_VARS_commID )
333 call atm_mesh%Create_communicator( &
335 this%QTRCVARS_manager, &
336 this%QTRC_VARS(1:qa), &
337 this%QTRC_VARS_commID )
340 call atm_mesh%Create_communicator( &
342 this%AUXVARS_manager, &
344 this%AUX_VARS_commID )
349 log_info(
"ATMOS_vars_setup",*)
'List of prognostic variables (ATMOS) '
350 log_info_cont(
'(1x,A,A24,A,A48,A,A12,A)') &
351 ' |',
'VARNAME ',
'|', &
352 'DESCRIPTION ',
'[',
'UNIT ',
']'
354 log_info_cont(
'(1x,A,I3,A,A24,A,A48,A,A12,A)') &
355 'NO.',iv,
'|',prgvar_info(iv)%NAME,
'|', prgvar_info(iv)%DESC,
'[', prgvar_info(iv)%UNIT,
']'
358 log_info_cont(
'(1x,A,I3,A,A24,A,A48,A,A12,A)') &
359 'NO.',
prgvar_num+iv,
'|',tracer_name(iv),
'|', tracer_desc(iv),
'[', tracer_unit(iv),
']'
365 call this%AUXVARS2D_manager%Init()
368 reg_file_hist = .true.
370 call this%AUXVARS2D_manager%Regist( &
372 this%AUX_VARS2D(iv), &
373 reg_file_hist, fill_zero=.true. )
380 call diagvar_manager%Init()
383 reg_file_hist = .true.
385 call diagvar_manager%Regist( &
387 diag_vars3d(iv), reg_file_hist )
389 this%DIAGVARS3D_HISTID(iv) = diag_vars3d(iv)%hist_id
391 call diagvar_manager%Final()
394 call diagvar_manager%Init()
397 reg_file_hist = .true.
399 call diagvar_manager%Regist( &
401 diag_vars2d(iv), reg_file_hist )
403 this%DIAGVARS2D_HISTID(iv) = diag_vars2d(iv)%hist_id
405 call diagvar_manager%Final()
409 is_specified = .true.
412 read(io_fid_conf,nml=param_atmos_vars_restart,iostat=ierr)
414 log_info(
"AtmosVars_Init",*)
'Not found namelist PARAM_ATMOS_VARS_RESTART. Default used.'
415 is_specified = .false.
416 elseif( ierr > 0 )
then
417 log_error(
"AtmosVars_Init",*)
'Not appropriate names in namelist PARAM_ATMOS_VARS_RESTART. Check!'
420 log_nml(param_atmos_vars_restart)
422 if (is_specified)
then
423 call atm_mesh%Setup_restartfile( this%restart_file, &
424 in_basename, in_postfix_timelabel, &
425 out_basename, out_postfix_timelabel, out_dtype, out_title, &
428 call atm_mesh%Setup_restartfile( this%restart_file, &
434 call monitor_reg(
'QTOT',
'water mass',
'kg', &
435 dv_monit_id(im_qtot), &
436 dim_type=
'ATM3D', is_tendency=.false. )
438 call monitor_reg(
'ENGT',
'total energy',
'J', &
439 dv_monit_id(im_engt), &
440 dim_type=
'ATM3D', is_tendency=.false. )
441 call monitor_reg(
'ENGP',
'potential energy',
'J', &
442 dv_monit_id(im_engp), &
443 dim_type=
'ATM3D', is_tendency=.false. )
444 call monitor_reg(
'ENGK',
'kinetic energy',
'J', &
445 dv_monit_id(im_engk), &
446 dim_type=
'ATM3D', is_tendency=.false. )
447 call monitor_reg(
'ENGI',
'internal energy',
'J', &
448 dv_monit_id(im_engi), &
449 dim_type=
'ATM3D', is_tendency=.false. )
454 this%check_range = check_range
455 this%check_total = check_total
456 log_info(
"ATMOS_vars_setup",*)
'Check value range of variables? : ', check_range
457 log_info(
"ATMOS_vars_setup",*)
'Check total value of variables? : ', check_total
460 end subroutine atmosvars_init
465 subroutine atmosvars_final( this )
471 log_info(
'AtmosVars_Final',*)
473 call this%restart_file%Final()
475 call this%PROGVARS_manager%Final()
476 deallocate( this%PROG_VARS )
478 call this%QTRCVARS_manager%Final()
479 deallocate( this%QTRC_VARS )
481 call this%AUXVARS_manager%Final()
482 deallocate( this%AUX_VARS )
484 call this%AUXVARS2D_manager%Final()
485 deallocate( this%AUX_VARS2D )
487 call this%PHYTENDS_manager%Final()
488 deallocate( this%PHY_TEND )
490 deallocate( this%DIAGVARS3D_HISTID )
493 end subroutine atmosvars_final
496 subroutine atmosvars_regist_physvar_manager( this, &
497 mp_AUXVARS2D_manager )
500 class(
atmosvars),
target,
intent(inout) :: this
504 this%ptr_MP_AUXVARS2D_manager => mp_auxvars2d_manager
507 end subroutine atmosvars_regist_physvar_manager
512 subroutine atmosvars_history( this )
527 mesh3d => this%PROG_VARS(1)%mesh
528 call mesh3d%GetMesh2D(mesh2d)
531 hst_id = this%PROG_VARS(v)%hist_id
536 hst_id = this%QTRC_VARS(v)%hist_id
540 call this%Calc_diagnostics()
543 hst_id = this%AUX_VARS(v)%hist_id
547 hst_id = this%AUX_VARS2D(v)%hist_id
551 do v = 1, this%PHYTEND_NUM_TOT
552 hst_id = this%PHY_TEND(v)%hist_id
559 call tmp_field3d%Init(
"tmp_field",
"", mesh3d)
561 hst_id = this%DIAGVARS3D_HISTID(v)
562 if ( hst_id > 0 )
then
567 call tmp_field3d%Final()
570 call tmp_field2d%Init(
"tmp_field",
"", mesh2d)
572 hst_id = this%DIAGVARS2D_HISTID(v)
573 if ( hst_id > 0 )
then
578 call tmp_field2d%Final()
581 end subroutine atmosvars_history
586 subroutine atmosvar_read_restart_file( this, atmos_mesh, dyncore )
596 class(
atmosvars),
intent(inout),
target :: this
597 class(
atmosmesh),
intent(in) :: atmos_mesh
610 log_info(
"ATMOSVar_read_restart_file",*)
'Open restart file (ATMOS) '
613 call this%restart_file%Open()
618 call this%restart_file%Read_var( dimtype_xyz, this%PROG_VARS(iv)%varname, &
622 call this%restart_file%Read_var( dimtype_xyz, this%AUX_VARS(iv)%varname, &
626 call this%restart_file%Read_var( dimtype_xyz, this%QTRC_VARS(iv)%varname, &
631 log_info(
"ATMOSVar_read_restart_file",*)
'Close restart file (ATMOS) '
632 call this%restart_file%Close()
637 call vars_calc_specific_heat( this )
640 call dyncore%update_therm_hyd( this%AUXVARS_manager )
644 this%PROGVARS_manager, this%AUXVARS_manager )
648 mesh3d => phyd_ref%mesh
649 do domid=1, mesh3d%LOCAL_MESH_NUM
650 lcmesh3d => mesh3d%lcmesh_list(domid)
652 do ke=lcmesh3d%NeS, lcmesh3d%NeE
653 phyd_ref%local(domid)%val(:,ke) = 0.0_rp
658 call this%Check( force = .true. )
661 call this%Calc_diagnostics()
664 call this%AUXVARS_manager%MeshFieldComm_Exchange()
668 mesh3d, atmos_mesh%element3D_operation )
671 end subroutine atmosvar_read_restart_file
676 subroutine atmosvar_write_restart_file( this )
686 integer :: iv, rf_vid
690 log_info(
"ATMOSVar_write_restart_file",*)
'Create restart file (ATMOS) '
693 call this%Check( force = .true. )
696 call this%restart_file%Create()
697 call prc_mpibarrier()
705 call this%restart_file%Def_var( this%PROG_VARS(iv), &
706 prgvar_info(iv)%DESC, rf_vid, dimtype_xyz )
710 call this%restart_file%Def_var( this%AUX_VARS(iv), &
711 auxvar_info(iv)%DESC, rf_vid, dimtype_xyz )
715 call this%restart_file%Def_var( this%QTRC_VARS(iv), &
716 tracer_desc(iv), rf_vid, dimtype_xyz )
719 call this%restart_file%End_def()
724 call this%restart_file%Write_var(rf_vid, this%PROG_VARS(iv) )
728 call this%restart_file%Write_var(rf_vid, this%AUX_VARS(iv) )
732 call this%restart_file%Write_var(rf_vid, this%QTRC_VARS(iv) )
736 log_info(
"ATMOSVar_write_restart_file",*)
'Close restart file (ATMOS) '
737 call this%restart_file%Close()
740 end subroutine atmosvar_write_restart_file
746 subroutine atmosvars_check( this, force )
754 logical,
intent(in),
optional :: force
765 character(len=H_MID) :: varname
771 if (
present(force) )
then
774 check = this%check_range
779 if ( iv == prgvar_therm_id ) cycle
781 mesh3d => this%PROG_VARS(iv)%mesh
782 do n=1, mesh3d%LOCAL_MESH_NUM
783 lcmesh => mesh3d%lcmesh_list(n)
784 elem => lcmesh%refElem
785 call this%PROG_VARS(iv)%GetLocalMeshField(n, lcfield)
786 write(varname,
'(a,i3.3,a)') this%PROG_VARS(iv)%varname//
'(domID=', n,
')'
787 call valcheck( elem%Np, 1, elem%Np, lcmesh%NeA, lcmesh%NeS, lcmesh%NeE, lcfield%val(:,:), &
792 mesh3d => this%PROG_VARS(1)%mesh
796 call atmosvars_calcdiagvar( this, vel_fields(iv)%varname, vel_fields(iv) )
800 call vel_fields(iv)%Final()
804 if (
present(force) )
then
807 check = this%check_total
810 mesh3d => this%PROG_VARS(1)%mesh
811 call work%Init(
"tmp",
"", mesh3d)
816 end subroutine atmosvars_check
821 subroutine atmosvars_monitor( this )
837 mesh3d => this%PROG_VARS(1)%mesh
838 call work%Init(
"tmp",
"", mesh3d)
845 if ( this%QTRC_VARS(iv)%monitor_id > 0 )
then
846 do n=1, mesh3d%LOCAL_MESH_NUM
847 lcmesh => mesh3d%lcmesh_list(n)
849 do ke=lcmesh%NeS, lcmesh%NeE
850 work%local(n)%val(:,ke) = ( this%AUX_VARS(auxvar_denshydro_id)%local(n)%val(:,ke) &
851 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) &
852 ) * this%QTRC_VARS(iv)%local(n)%val(:,ke)
858 if ( dv_monit_id(im_qtot) > 0 )
then
859 do n=1, mesh3d%LOCAL_MESH_NUM
860 lcmesh => mesh3d%lcmesh_list(n)
862 do ke=lcmesh%NeS, lcmesh%NeE
863 work%local(n)%val(:,ke) = ( this%AUX_VARS(auxvar_denshydro_id)%local(n)%val(:,ke) &
864 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) &
865 ) * ( 1.0_rp - this%AUX_VARS(auxvar_qdry_id)%local(n)%val(:,ke) )
873 if ( dv_monit_id(im_engt) > 0 )
then
874 call atmosvars_calcdiagvar( this,
'ENGT', work )
877 if ( dv_monit_id(im_engp) > 0 )
then
878 call atmosvars_calcdiagvar( this,
'ENGP', work )
881 if ( dv_monit_id(im_engk) > 0 )
then
882 call atmosvars_calcdiagvar( this,
'ENGK', work )
885 if ( dv_monit_id(im_engi) > 0 )
then
886 call atmosvars_calcdiagvar( this,
'ENGI', work )
893 end subroutine atmosvars_monitor
900 var, DENS_hyd, PRES_hyd, lcmesh3D )
904 integer,
intent(in) :: domid
908 integer,
intent(in) :: varid
911 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
918 call prgvars_list%Get(varid, field)
919 call field%GetLocalMeshField(domid, var)
921 if (
present(dens_hyd))
then
922 call auxvars_list%Get(auxvar_denshydro_id, field)
923 call field%GetLocalMeshField(domid, dens_hyd)
925 if (
present(pres_hyd))
then
926 call auxvars_list%Get(auxvar_preshydro_id, field)
927 call field%GetLocalMeshField(domid, pres_hyd)
930 if (
present(lcmesh3d))
then
931 call mesh%GetLocalMesh( domid, lcmesh )
936 if (
present(lcmesh3d)) lcmesh3d => lcmesh
945 DDENS, MOMX, MOMY, MOMZ, THERM, &
946 DENS_hyd, PRES_hyd, Rtot, CVtot, CPtot, &
949 integer,
intent(in) :: domid
956 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
963 call prgvars_list%Get(prgvar_ddens_id, field)
964 call field%GetLocalMeshField(domid, ddens)
966 call prgvars_list%Get(prgvar_momx_id, field)
967 call field%GetLocalMeshField(domid, momx)
969 call prgvars_list%Get(prgvar_momy_id, field)
970 call field%GetLocalMeshField(domid, momy)
972 call prgvars_list%Get(prgvar_momz_id, field)
973 call field%GetLocalMeshField(domid, momz)
975 call prgvars_list%Get(prgvar_therm_id, field)
976 call field%GetLocalMeshField(domid, therm)
979 call auxvars_list%Get(auxvar_denshydro_id, field)
980 call field%GetLocalMeshField(domid, dens_hyd)
982 call auxvars_list%Get(auxvar_preshydro_id, field)
983 call field%GetLocalMeshField(domid, pres_hyd)
985 call auxvars_list%Get(auxvar_rtot_id, field)
986 call field%GetLocalMeshField(domid, rtot)
988 call auxvars_list%Get(auxvar_cvtot_id, field)
989 call field%GetLocalMeshField(domid, cvtot)
991 call auxvars_list%Get(auxvar_cptot_id, field)
992 call field%GetLocalMeshField(domid, cptot)
996 if (
present(lcmesh3d) )
then
997 call mesh%GetLocalMesh( domid, lcmesh )
1002 if (
present(lcmesh3d)) lcmesh3d => lcmesh
1011 PREC, PREC_ENGI, lcmesh2D )
1014 integer,
intent(in) :: domid
1015 class(
meshbase),
intent(in) :: mesh
1018 class(
localmesh2d),
pointer,
intent(out),
optional :: lcmesh2d
1026 call field%GetLocalMeshField(domid, prec)
1029 call field%GetLocalMeshField(domid, prec_engi)
1031 if (
present(lcmesh2d))
then
1032 call mesh%GetLocalMesh( domid, lcmesh )
1037 if (
present(lcmesh2d)) lcmesh2d => lcmesh
1050 integer,
intent(in) :: domid
1051 class(
meshbase),
intent(in) :: mesh
1053 integer,
intent(in) :: varid
1055 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
1062 call trcvars_list%Get(varid, field)
1063 call field%GetLocalMeshField(domid, var)
1065 if (
present(lcmesh3d))
then
1066 call mesh%GetLocalMesh( domid, lcmesh )
1071 if (
present(lcmesh3d)) lcmesh3d => lcmesh
1080 var, var_tp, lcmesh3D )
1082 use scale_atmos_hydrometeor,
only: &
1083 atmos_hydrometeor_dry, &
1086 integer,
intent(in) :: domid
1087 class(
meshbase),
intent(in) :: mesh
1092 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
1097 integer :: iq, tend_iq
1102 if ( atmos_hydrometeor_dry )
then
1103 iq = 0; tend_iq = phytend_num1+1
1105 iq = i_qv; tend_iq = phytend_num1 + i_qv
1108 call trcvars_list%Get(iq, field)
1109 call field%GetLocalMeshField(domid, var)
1111 call forcing_list%Get(tend_iq, field)
1112 call field%GetLocalMeshField(domid, var_tp)
1114 if (
present(lcmesh3d))
then
1115 call mesh%GetLocalMesh( domid, lcmesh )
1120 if (
present(lcmesh3d)) lcmesh3d => lcmesh
1130 var_list, lcmesh3D )
1133 integer,
intent(in) :: domid
1134 class(
meshbase),
intent(in) :: mesh
1136 integer,
intent(in) :: varid_s
1138 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
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)
1151 if (
present(lcmesh3d))
then
1152 call mesh%GetLocalMesh( domid, lcmesh )
1157 if (
present(lcmesh3d)) lcmesh3d => lcmesh
1170 integer,
intent(in) :: domid
1171 class(
meshbase),
intent(in) :: mesh
1174 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
1181 call phyauxvars_list%Get(auxvar_pres_id, field)
1182 call field%GetLocalMeshField(domid, pres)
1184 call phyauxvars_list%Get(auxvar_pt_id, field)
1185 call field%GetLocalMeshField(domid, pt)
1189 if (
present(lcmesh3d))
then
1190 call mesh%GetLocalMesh( domid, lcmesh )
1195 if (
present(lcmesh3d)) lcmesh3d => lcmesh
1204 DENS_tp, MOMX_tp, MOMY_tp, MOMZ_tp, RHOT_tp, RHOH_p, &
1209 integer,
intent(in) :: domid
1210 class(
meshbase),
intent(in) :: mesh
1212 class(
localmeshfieldbase),
pointer,
intent(out) :: dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp
1215 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
1224 call phytends_list%Get(phytend_dens_id, field)
1225 call field%GetLocalMeshField(domid, dens_tp)
1227 call phytends_list%Get(phytend_momx_id, field)
1228 call field%GetLocalMeshField(domid, momx_tp)
1230 call phytends_list%Get(phytend_momy_id, field)
1231 call field%GetLocalMeshField(domid, momy_tp)
1233 call phytends_list%Get(phytend_momz_id, field)
1234 call field%GetLocalMeshField(domid, momz_tp)
1236 call phytends_list%Get(phytend_rhot_id, field)
1237 call field%GetLocalMeshField(domid, rhot_tp)
1239 call phytends_list%Get(phytend_rhoh_id, field)
1240 call field%GetLocalMeshField(domid, rhoh_p)
1242 if (
present(rhoq_tp) )
then
1244 call phytends_list%Get(phytend_num1+iq, field)
1245 call field%GetLocalMeshField(domid, rhoq_tp(iq)%ptr)
1250 if (
present(lcmesh3d) )
then
1251 call mesh%GetLocalMesh( domid, lcmesh )
1256 if (
present(lcmesh3d) ) lcmesh3d => lcmesh
1269 integer,
intent(in) :: domid
1270 class(
meshbase),
intent(in) :: mesh
1272 integer,
intent(in) :: qtrcid
1279 call phytends_list%Get(phytend_num1 + qtrcid, field)
1280 call field%GetLocalMeshField(domid, rhoq_tp)
1288 subroutine atmosvars_calculatediagnostics( this )
1289 use scale_const,
only: &
1290 rdry => const_rdry, &
1291 cpdry => const_cpdry, &
1292 cvdry => const_cvdry, &
1293 pres00 => const_pre00
1294 use scale_tracer,
only: &
1295 tracer_mass, tracer_r, tracer_cv, tracer_cp
1296 use scale_atmos_thermodyn,
only: &
1297 atmos_thermodyn_specific_heat
1299 class(
atmosvars),
intent(inout),
target :: this
1313 call vars_calc_specific_heat( this )
1316 do varid=auxvar_thermhydro_id+1, auxvar_pt_id
1317 field => this%AUX_VARS(varid)
1318 do n=1, field%mesh%LOCAL_MESH_NUM
1320 field%mesh, this%QTRCVARS_manager, &
1323 elem3d => lcmesh3d%refElem3D
1325 call vars_calc_diagnosevar_lc( &
1326 field%varname, field%local(n)%val, &
1327 this%PROG_VARS(prgvar_ddens_id)%local(n)%val, &
1328 this%PROG_VARS(prgvar_momx_id)%local(n)%val, &
1329 this%PROG_VARS(prgvar_momy_id)%local(n)%val, &
1330 this%PROG_VARS(prgvar_momz_id)%local(n)%val, &
1331 this%AUX_VARS(auxvar_pres_id)%local(n)%val, &
1332 this%AUX_VARS(auxvar_qdry_id)%local(n)%val, &
1334 this%AUX_VARS(auxvar_denshydro_id)%local(n)%val, &
1335 this%AUX_VARS(auxvar_preshydro_id)%local(n)%val, &
1336 this%AUX_VARS(auxvar_rtot_id )%local(n)%val, &
1337 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val, &
1338 this%AUX_VARS(auxvar_cptot_id)%local(n)%val, &
1339 lcmesh3d, lcmesh3d%refElem3D )
1344 end subroutine atmosvars_calculatediagnostics
1347 subroutine atmosvars_calcdiagvar( this, field_name, field_work )
1348 use scale_const,
only: &
1349 rdry => const_rdry, &
1350 cpdry => const_cpdry, &
1351 cvdry => const_cvdry, &
1352 pres00 => const_pre00
1353 use scale_tracer,
only: &
1354 tracer_mass, tracer_r, tracer_cv, tracer_cp
1355 use scale_atmos_thermodyn,
only: &
1356 atmos_thermodyn_specific_heat
1360 character(*),
intent(in) :: field_name
1377 if ( field_name ==
'Umet' )
then
1378 is_uvmet = .true.; uvmet_i = 1
1379 else if ( field_name ==
'Vmet' )
then
1380 is_uvmet = .true.; uvmet_i = 2
1383 field_work%varname = field_name
1385 do n=1, field_work%mesh%LOCAL_MESH_NUM
1387 field_work%mesh, this%QTRCVARS_manager, &
1390 if ( .not. is_uvmet )
then
1391 call vars_calc_diagnosevar_lc( field_name, field_work%local(n)%val, &
1392 this%PROG_VARS(prgvar_ddens_id)%local(n)%val, &
1393 this%PROG_VARS(prgvar_momx_id)%local(n)%val, &
1394 this%PROG_VARS(prgvar_momy_id)%local(n)%val, &
1395 this%PROG_VARS(prgvar_momz_id)%local(n)%val, &
1396 this%AUX_VARS(auxvar_pres_id)%local(n)%val, &
1397 this%AUX_VARS(auxvar_qdry_id)%local(n)%val, &
1399 this%AUX_VARS(auxvar_denshydro_id)%local(n)%val, &
1400 this%AUX_VARS(auxvar_preshydro_id)%local(n)%val, &
1401 this%AUX_VARS(auxvar_rtot_id )%local(n)%val, &
1402 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val, &
1403 this%AUX_VARS(auxvar_cptot_id)%local(n)%val, &
1404 lcmesh3d, lcmesh3d%refElem3D )
1406 call field_work_uvmet(1)%Init(
'Umet',
'', field_work%mesh )
1407 call field_work_uvmet(2)%Init(
'Vmet',
'', field_work%mesh )
1408 call this%mesh%Calc_UVmet( &
1409 this%PROG_VARS(prgvar_momx_id), this%PROG_VARS(prgvar_momy_id), &
1410 field_work_uvmet(1), field_work_uvmet(2) )
1412 do ke=lcmesh3d%NeS, lcmesh3d%NeE
1413 field_work%local(n)%val(:,ke) = field_work_uvmet(uvmet_i)%local(n)%val(:,ke) &
1414 / ( this%AUX_VARS (auxvar_denshydro_id)%local(n)%val(:,ke) &
1415 + this%PROG_VARS(prgvar_ddens_id )%local(n)%val(:,ke) )
1417 call field_work_uvmet(1)%Final()
1418 call field_work_uvmet(2)%Final()
1423 end subroutine atmosvars_calcdiagvar
1426 subroutine atmosvars_calcdiagvar2d( this, field_name, field_work )
1429 character(*),
intent(in) :: field_name
1430 type(
meshfield2d),
intent(inout),
target :: field_work
1438 field_work%varname = field_name
1440 do n=1, field_work%mesh%LOCAL_MESH_NUM
1441 lcmesh2d => field_work%mesh%lcmesh_list(n)
1442 call vars_calc_diagnosevar2d_lc( field_name, field_work%local(n)%val, &
1443 this%ptr_MP_AUXVARS2D_manager, &
1444 field_work%mesh, lcmesh2d, lcmesh2d%refElem2D )
1448 end subroutine atmosvars_calcdiagvar2d
1453 subroutine vars_calc_diagnosevar_lc( field_name, var_out, &
1454 DDENS_, MOMX_, MOMY_, MOMZ_, PRES_, QDRY_, QTRC, &
1455 DENS_hyd, PRES_hyd, Rtot, CVtot, CPTot, &
1458 use scale_const,
only: &
1459 grav => const_grav, &
1460 rdry => const_rdry, &
1461 rvap => const_rvap, &
1462 cpdry => const_cpdry, &
1463 cvdry => const_cvdry, &
1464 pres00 => const_pre00
1465 use scale_tracer,
only: &
1467 tracer_cv, tracer_engi0
1468 use scale_atmos_saturation,
only: &
1469 atmos_saturation_psat_liq
1474 character(*),
intent(in) :: field_name
1475 real(rp),
intent(out) :: var_out(elem%np,lcmesh%nea)
1476 real(rp),
intent(in) :: ddens_(elem%np,lcmesh%nea)
1477 real(rp),
intent(in) :: momx_(elem%np,lcmesh%nea)
1478 real(rp),
intent(in) :: momy_(elem%np,lcmesh%nea)
1479 real(rp),
intent(in) :: momz_(elem%np,lcmesh%nea)
1480 real(rp),
intent(in) :: pres_(elem%np,lcmesh%nea)
1481 real(rp),
intent(in) :: qdry_(elem%np,lcmesh%nea)
1483 real(rp),
intent(in) :: dens_hyd(elem%np,lcmesh%nea)
1484 real(rp),
intent(in) :: pres_hyd(elem%np,lcmesh%nea)
1485 real(rp),
intent(in) :: rtot (elem%np,lcmesh%nea)
1486 real(rp),
intent(in) :: cvtot(elem%np,lcmesh%nea)
1487 real(rp),
intent(in) :: cptot(elem%np,lcmesh%nea)
1491 real(rp) :: dens(elem%np), temp(elem%np)
1492 real(rp) :: mom_u1(elem%np), mom_u2(elem%np), g_11(elem%np), g_12(elem%np), g_22(elem%np)
1493 real(rp) :: psat(elem%np)
1498 select case(trim(field_name))
1502 var_out(:,ke) = ddens_(:,ke) + dens_hyd(:,ke)
1508 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1509 var_out(:,ke) = momx_(:,ke) / dens(:)
1515 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1516 var_out(:,ke) = momy_(:,ke) / dens(:)
1522 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1523 var_out(:,ke) = momz_(:,ke) / dens(:)
1530 var_out(:,ke) = pres_(:,ke) - pres_hyd(:,ke)
1536 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) )
1542 var_out(:,ke) = pres_(:,ke) / ( rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) ) &
1543 - pres_hyd(:,ke) / ( rdry * dens_hyd(:,ke) )
1549 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1550 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * dens(:) ) * ( pres00 / pres_(:,ke) )**( rtot(:,ke) / cptot(:,ke) )
1556 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1557 var_out(:,ke) = pres_(:,ke) / (rtot(:,ke) * dens(:) ) * ( pres00 / pres_(:,ke) )**( rtot(:,ke) / cptot(:,ke) ) &
1558 - pres00/rdry * (pres_hyd(:,ke)/pres00)**(cvdry/cpdry) / dens_hyd(:,ke)
1562 if ( atmos_hydrometeor_dry )
then
1563 var_out(:,ke) = 0.0_rp
1565 call tracer_inq_id(
"QV", iq_qv )
1569 temp(:) = pres_(:,ke) / (rtot(:,ke) * (ddens_(:,ke) + dens_hyd(:,ke)) )
1571 call atmos_saturation_psat_liq( &
1572 elem%Np, 1, elem%Np, temp(:), &
1575 var_out(:,ke) = ( ddens_(:,ke) + dens_hyd(:,ke) ) * qtrc(iq_qv)%ptr%val(:,ke) &
1576 / psat(:) * rvap * temp(:) * 100.0_rp
1583 ke2d = lcmesh%EMap3Dto2D(ke)
1585 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1586 g_11(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,1)
1587 g_12(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,2)
1588 g_22(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,2)
1590 mom_u1(:) = g_11(:) * momx_(:,ke) + g_12(:) * momy_(:,ke)
1591 mom_u2(:) = g_12(:) * momx_(:,ke) + g_22(:) * momy_(:,ke)
1593 var_out(:,ke) = 0.5_rp * ( momx_(:,ke) * mom_u1(:) + momy_(:,ke) * mom_u2(:) + momz_(:,ke)**2 ) / dens(:)
1599 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1600 var_out(:,ke) = dens(:) * grav * lcmesh%zlev(:,ke)
1606 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1607 var_out(:,ke) = qdry_(:,ke) * pres_(:,ke) / rtot(:,ke) * cvdry
1609 var_out(:,ke) = var_out(:,ke) &
1610 + qtrc(iq)%ptr%val(:,ke) * ( pres_(:,ke) / rtot(:,ke) * tracer_cv(iq) + dens(:) * tracer_engi0(iq) )
1617 ke2d = lcmesh%EMap3Dto2D(ke)
1619 dens(:) = ddens_(:,ke) + dens_hyd(:,ke)
1621 g_11(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,1)
1622 g_12(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,1,2)
1623 g_22(:) = lcmesh%G_ij(elem%IndexH2Dto3D,ke2d,2,2)
1624 mom_u1(:) = g_11(:) * momx_(:,ke) + g_12(:) * momy_(:,ke)
1625 mom_u2(:) = g_12(:) * momx_(:,ke) + g_22(:) * momy_(:,ke)
1628 var_out(:,ke) = qdry_(:,ke) * pres_(:,ke) / rtot(:,ke) * cvdry
1630 var_out(:,ke) = var_out(:,ke) &
1631 + qtrc(iq)%ptr%val(:,ke) * ( pres_(:,ke) / rtot(:,ke) * tracer_cv(iq) + dens(:) * tracer_engi0(iq) )
1635 0.5_rp * ( momx_(:,ke) * mom_u1(:) + momy_(:,ke) * mom_u2(:) + momz_(:,ke)**2 ) / dens(:) &
1637 + dens(:) * grav * lcmesh%pos_en(:,ke,3)
1641 log_error(
"AtmosVars_calc_diagnoseVar_lc",*)
'The name of diagnostic variable is not suported. Check!', field_name
1647 end subroutine vars_calc_diagnosevar_lc
1650 subroutine vars_calc_diagnosevar2d_lc( field_name, & ! (in)
1652 mp_auxvars2d, mesh2d, lcmesh, elem )
1660 character(*),
intent(in) :: field_name
1661 real(rp),
intent(out) :: var_out(elem%np,lcmesh%nea)
1670 select case(trim(field_name))
1671 case(
'RAIN',
'SNOW')
1673 lcmesh%lcdomID, mesh2d, mp_auxvars2d, &
1674 sflx_rain_mp, sflx_snow_mp, sflx_engi_mp )
1677 select case(trim(field_name))
1680 do ke=lcmesh%NeS, lcmesh%NeE
1681 var_out(:,ke) = sflx_rain_mp%val(:,ke)
1685 do ke=lcmesh%NeS, lcmesh%NeE
1686 var_out(:,ke) = sflx_snow_mp%val(:,ke)
1689 log_error(
"AtmosVars_calc_diagnoseVar2D_lc",*)
'The name of diagnostic variable is not suported. Check!', field_name
1694 end subroutine vars_calc_diagnosevar2d_lc
1697 subroutine vars_calc_specific_heat( this )
1698 use scale_const,
only: &
1699 rdry => const_rdry, &
1700 cpdry => const_cpdry, &
1701 cvdry => const_cvdry, &
1702 pres00 => const_pre00
1703 use scale_tracer,
only: &
1704 tracer_mass, tracer_r, tracer_cv, tracer_cp
1705 use scale_atmos_thermodyn,
only: &
1706 atmos_thermodyn_specific_heat
1708 class(
atmosvars),
intent(inout),
target :: this
1718 real(rp),
allocatable :: q_tmp(:,:)
1722 do n=1, this%AUX_VARS(1)%mesh%LOCAL_MESH_NUM
1723 lcmesh3d => this%AUX_VARS(1)%mesh%lcmesh_list(n)
1724 elem3d => lcmesh3d%refElem3D
1725 allocate( q_tmp(elem3d%Np,qa) )
1728 do ke = lcmesh3d%NeS, lcmesh3d%NeE
1730 q_tmp(:,iq) = this%QTRC_VARS(iq)%local(n)%val(:,ke)
1732 call atmos_thermodyn_specific_heat( &
1733 elem3d%Np, 1, elem3d%Np, qa, &
1734 q_tmp(:,:), tracer_mass(:), tracer_r(:), tracer_cv(:), tracer_cp(:), &
1735 this%AUX_VARS(auxvar_qdry_id )%local(n)%val(:,ke), &
1736 this%AUX_VARS(auxvar_rtot_id )%local(n)%val(:,ke), &
1737 this%AUX_VARS(auxvar_cvtot_id)%local(n)%val(:,ke), &
1738 this%AUX_VARS(auxvar_cptot_id)%local(n)%val(:,ke) )
1742 end subroutine vars_calc_specific_heat
module Atmosphere / Physics / Cloud Microphysics
subroutine, public atmosphympvars_getlocalmeshfields_sfcflx(domid, mesh, sfcflx_list, sflx_rain, sflx_snow, sflx_engi)
module Atmosphere / Physics / Variables
integer, parameter, public atmos_diagvars2d_num
real(rp), dimension(prgvar_num), parameter progvars_check_min
integer, parameter, public atmos_diagvars_dens_id
integer, parameter, public atmos_diagvars_w_id
subroutine, public atmosvars_getlocalmeshsfcvar(domid, mesh, auxvars2d_list, prec, prec_engi, lcmesh2d)
type(variableinfo), dimension(atmos_diagvars3d_num), public atmos_diagvars3d_vinfo
subroutine, public atmosvars_getlocalmeshqtrc_qv(domid, mesh, trcvars_list, forcing_list, var, var_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshprgvar(domid, mesh, prgvars_list, auxvars_list, varid, var, dens_hyd, pres_hyd, lcmesh3d)
integer, parameter, public atmos_auxvars2d_prec_id
integer, parameter, public atmos_diagvars3d_num
subroutine, public atmosvars_getlocalmeshprgvars(domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
integer, parameter, public atmos_diagvars_vmet_id
integer, parameter, public atmos_diagvars_rh_id
integer, parameter, public atmos_diagvars_engi_id
integer, parameter, public atmos_diagvars_umet_id
integer, parameter, public atmos_diagvars_qdry_id
subroutine, public atmosvars_getlocalmeshqtrcvarlist(domid, mesh, trcvars_list, varid_s, var_list, lcmesh3d)
integer, parameter, public atmos_auxvars2d_prec_engi_id
integer, parameter, public atmos_diagvars_engk_id
type(variableinfo), dimension(atmos_diagvars2d_num), public atmos_diagvars2d_vinfo
integer, parameter, public atmos_diagvars_v_id
integer, parameter, public atmos_auxvars2d_num
type(variableinfo), dimension(atmos_auxvars2d_num), public atmos_auxvars2d_vinfo
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)
integer, parameter, public atmos_diagvars_t_id
subroutine, public atmosvars_getlocalmeshqtrcvar(domid, mesh, trcvars_list, varid, var, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphyauxvars(domid, mesh, phyauxvars_list, pres, pt, lcmesh3d)
integer, parameter, public atmos_diagvars_snow_id
integer, parameter, public atmos_diagvars_engt_id
integer, parameter, public atmos_diagvars_engp_id
integer, parameter, public atmos_diagvars_rain_id
integer, parameter, public atmos_diagvars_u_id
module FElib / Fluid dyn solver / Atmosphere / driver (3D nonhydrostatic model)
module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common
integer, parameter, public prgvar_scalar_num
integer, parameter, public phytend_momz_id
integer, parameter, public prgvar_momy_id
integer, parameter, public auxvar_cptot_id
integer, parameter, public prgvar_therm_id
integer, parameter, public auxvar_qdry_id
integer, parameter, public auxvar_cvtot_id
integer, parameter, public prgvar_ddens_id
integer, parameter, public prgvar_momz_id
integer, parameter, public phytend_rhot_id
integer, parameter, public phytend_momx_id
integer, parameter, public phytend_num
integer, parameter, public phytend_momy_id
integer, parameter, public prgvar_momx_id
integer, parameter, public auxvar_num
subroutine, public atm_dyn_dgm_nonhydro3d_common_get_varinfo(prgvar_info, auxvar_info, phytend_info)
integer, parameter, public auxvar_pt_id
integer, parameter, public prgvar_hvec_num
integer, parameter, public prgvar_num
integer, parameter, public phytend_rhoh_id
integer, parameter, public auxvar_preshydro_id
integer, parameter, public auxvar_thermhydro_id
integer, parameter, public auxvar_denshydro_id
integer, parameter, public phytend_dens_id
integer, parameter, public auxvar_rtot_id
subroutine, public atm_dyn_dgm_nonhydro3d_common_setup_variables(prgvars, qtrcvars, auxvars, phytends, prgvar_manager, qtrcvar_manager, auxvar_manager, phytend_manager, phytend_num_tot, mesh3d, prgvar_varinfo)
integer, parameter, public auxvar_preshydro_ref_id
integer, parameter, public auxvar_pres_id
module FElib / Element / Base
module FElib / File / History
module FElib / File / Monitor
module FElib / File / Restart
type(file_restart_meshfield), public restart_file
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Data / base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
integer, public meshbase3d_dimtypeid_xyz
module FElib / Mesh / Base
module FElib / Data / base
module FElib / Data / Statistics
module FElib / Data / Communication base
module FElib / Data / Communication 3D cubic domain
FElib / model framework / variable manager.
Derived type to manage a computational mesh (base class)
Derived type to manage variables with atmospheric component.
Derived type to provide a driver of dynamical core with the atmospheric nonhydrostatic equations.
Derived type representing a 2D reference element.
Derived type representing a 3D reference element.
Derived type representing an arbitrary finite element.
Derived type to manage a local 3D computational domain.
Derived type to manage a local computational domain (base type)
Derived type representing a field with local mesh (base type)
Derived type representing a field with 2D mesh.
Derived type representing a field with 3D mesh.
Derived type representing a field (base type)
Container to save a pointer of MeshField(1D, 2D, 3D) object.
Base derived type to manage data communication with 3D cubic domain.