Skip to content

Commit

Permalink
Synced with NCAR repo.
Browse files Browse the repository at this point in the history
  • Loading branch information
dustinswales committed Jun 27, 2019
1 parent fa05574 commit f7915b9
Show file tree
Hide file tree
Showing 42 changed files with 1,273 additions and 1,084 deletions.
8 changes: 2 additions & 6 deletions physics/GFS_DCNV_generic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -139,8 +139,6 @@ end subroutine GFS_DCNV_generic_post_finalize
!! | ncnvcld3d | number_of_convective_3d_cloud_fields | number of convective 3d clouds fields | count | 0 | integer | | in | F |
!! | rainc | lwe_thickness_of_convective_precipitation_amount_on_dynamics_timestep | convective rain at this time step | m | 1 | real | kind_phys | inout | F |
!! | cldwrk | cumulative_cloud_work_function | cumulative cloud work function (valid only with sas) | m2 s-1 | 1 | real | kind_phys | inout | F |
!! | cnvprcp | cumulative_lwe_thickness_of_convective_precipitation_amount | cumulative convective precipitation | m | 1 | real | kind_phys | inout | F |
!! | cnvprcpb | cumulative_lwe_thickness_of_convective_precipitation_amount_in_bucket | cumulative convective precipitation in bucket | m | 1 | real | kind_phys | inout | F |
!! | dt3dt | cumulative_change_in_temperature_due_to_deep_convection | cumulative change in temperature due to deep conv. | K | 2 | real | kind_phys | inout | F |
!! | dq3dt | cumulative_change_in_water_vapor_specific_humidity_due_to_deep_convection | cumulative change in water vapor specific humidity due to deep conv. | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | du3dt | cumulative_change_in_x_wind_due_to_deep_convection | cumulative change in x wind due to deep convection | m s-1 | 2 | real | kind_phys | inout | F |
Expand Down Expand Up @@ -169,7 +167,7 @@ end subroutine GFS_DCNV_generic_post_finalize
subroutine GFS_DCNV_generic_post_run (im, levs, lssav, ldiag3d, lgocart, ras, cscnv, do_ca, &
isppt_deep, frain, rain1, dtf, cld1d, save_u, save_v, save_t, save_qv, gu0, gv0, gt0, &
gq0_water_vapor, ud_mf, dd_mf, dt_mf, con_g, clw_ice, clw_liquid, npdf3d, num_p3d, ncnvcld3d, &
rainc, cldwrk, cnvprcp, cnvprcpb, dt3dt, dq3dt, du3dt, dv3dt, upd_mf, dwn_mf, det_mf, dqdti, &
rainc, cldwrk, dt3dt, dq3dt, du3dt, dv3dt, upd_mf, dwn_mf, det_mf, dqdti, &
cnvqci, upd_mfi, dwn_mfi, det_mfi, cnvw, cnvc, cnvw_phy_f3d, cnvc_phy_f3d, &
cape, tconvtend, qconvtend, uconvtend, vconvtend, errmsg, errflg)

Expand All @@ -189,7 +187,7 @@ subroutine GFS_DCNV_generic_post_run (im, levs, lssav, ldiag3d, lgocart, ras, cs
real(kind=kind_phys), dimension(im,levs), intent(in) :: clw_ice, clw_liquid
integer, intent(in) :: npdf3d, num_p3d, ncnvcld3d

real(kind=kind_phys), dimension(im), intent(inout) :: rainc, cldwrk, cnvprcp, cnvprcpb
real(kind=kind_phys), dimension(im), intent(inout) :: rainc, cldwrk
! dt3dt, dq3dt, du3dt, dv3dt upd_mf, dwn_mf, det_mf only allocated if ldiag3d == .true.
real(kind=kind_phys), dimension(:,:), intent(inout) :: dt3dt, dq3dt, du3dt, dv3dt
real(kind=kind_phys), dimension(:,:), intent(inout) :: upd_mf, dwn_mf, det_mf
Expand Down Expand Up @@ -246,8 +244,6 @@ subroutine GFS_DCNV_generic_post_run (im, levs, lssav, ldiag3d, lgocart, ras, cs
if (lssav) then
do i=1,im
cldwrk (i) = cldwrk (i) + cld1d(i) * dtf
cnvprcp(i) = cnvprcp(i) + rainc(i)
cnvprcpb(i) = cnvprcpb(i) + rainc(i)
enddo

if (ldiag3d) then
Expand Down
87 changes: 53 additions & 34 deletions physics/GFS_MP_generic.F90
Original file line number Diff line number Diff line change
Expand Up @@ -104,6 +104,7 @@ end subroutine GFS_MP_generic_post_init
!! | imp_physics | flag_for_microphysics_scheme | choice of microphysics scheme | flag | 0 | integer | | in | F |
!! | imp_physics_gfdl | flag_for_gfdl_microphysics_scheme | choice of GFDL microphysics scheme | flag | 0 | integer | | in | F |
!! | imp_physics_thompson | flag_for_thompson_microphysics_scheme | choice of Thompson microphysics scheme | flag | 0 | integer | | in | F |
!! | imp_physics_mg | flag_for_morrison_gettelman_microphysics_scheme | choice of Morrison-Gettelman microphysics scheme | flag | 0 | integer | | in | F |
!! | cal_pre | flag_for_precipitation_type_algorithm | flag controls precip type algorithm | flag | 0 | logical | | in | F |
!! | lssav | flag_diagnostics | logical flag for storing diagnostics | flag | 0 | logical | | in | F |
!! | ldiag3d | flag_diagnostics_3D | flag for 3d diagnostic fields | flag | 0 | logical | | in | F |
Expand Down Expand Up @@ -140,10 +141,13 @@ end subroutine GFS_MP_generic_post_init
!! | doms_diag | dominant_snow_type | dominant snow type | none | 1 | real | kind_phys | inout | F |
!! | tprcp | nonnegative_lwe_thickness_of_precipitation_amount_on_dynamics_timestep | total precipitation amount in each time step | m | 1 | real | kind_phys | inout | F |
!! | srflag | flag_for_precipitation_type | snow/rain flag for precipitation | flag | 1 | real | kind_phys | inout | F |
!! | sr | ratio_of_snowfall_to_rainfall | snow ratio: ratio of snow to total precipitation | frac | 1 | real | kind_phys | in | F |
!! | cnvprcp | cumulative_lwe_thickness_of_convective_precipitation_amount | cumulative convective precipitation | m | 1 | real | kind_phys | inout | F |
!! | totprcp | accumulated_lwe_thickness_of_precipitation_amount | accumulated total precipitation | m | 1 | real | kind_phys | inout | F |
!! | totice | accumulated_lwe_thickness_of_ice_amount | accumulated ice precipitation | kg m-2 | 1 | real | kind_phys | inout | F |
!! | totsnw | accumulated_lwe_thickness_of_snow_amount | accumulated snow precipitation | kg m-2 | 1 | real | kind_phys | inout | F |
!! | totgrp | accumulated_lwe_thickness_of_graupel_amount | accumulated graupel precipitation | kg m-2 | 1 | real | kind_phys | inout | F |
!! | cnvprcpb | cumulative_lwe_thickness_of_convective_precipitation_amount_in_bucket | cumulative convective precipitation in bucket | m | 1 | real | kind_phys | inout | F |
!! | totprcpb | accumulated_lwe_thickness_of_precipitation_amount_in_bucket | accumulated total precipitation in bucket | m | 1 | real | kind_phys | inout | F |
!! | toticeb | accumulated_lwe_thickness_of_ice_amount_in_bucket | accumulated ice precipitation in bucket | kg m-2 | 1 | real | kind_phys | inout | F |
!! | totsnwb | accumulated_lwe_thickness_of_snow_amount_in_bucket | accumulated snow precipitation in bucket | kg m-2 | 1 | real | kind_phys | inout | F |
Expand Down Expand Up @@ -173,18 +177,19 @@ end subroutine GFS_MP_generic_post_init
!> \section gfs_mp_gen GFS MP Generic Post General Algorithm
!> @{
subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, ntrac, imp_physics, imp_physics_gfdl, &
imp_physics_thompson, cal_pre, lssav, ldiag3d, cplflx, cplchm, con_g, dtf, frain, rainc, rain1, rann, xlat, xlon, &
gt0, gq0, prsl, prsi, phii, tsfc, ice, snow, graupel, save_t, save_qv, rain0, ice0, snow0, graupel0, del, &
rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, srflag, totprcp, totice, totsnw, &
totgrp, totprcpb, toticeb, totsnwb, totgrpb, dt3dt, dq3dt, rain_cpl, rainc_cpl, snow_cpl, pwat, &
imp_physics_thompson, imp_physics_mg, cal_pre, lssav, ldiag3d, cplflx, cplchm, con_g, dtf, frain, rainc, rain1, &
rann, xlat, xlon, gt0, gq0, prsl, prsi, phii, tsfc, ice, snow, graupel, save_t, save_qv, rain0, ice0, snow0, &
graupel0, del, rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, srflag, sr, cnvprcp, totprcp, totice, &
totsnw, totgrp, cnvprcpb, totprcpb, toticeb, totsnwb, totgrpb, dt3dt, dq3dt, rain_cpl, rainc_cpl, snow_cpl, pwat, &
do_sppt, dtdtr, dtdtc, drain_cpl, dsnow_cpl, lsm, lsm_ruc, raincprv, rainncprv, iceprv, snowprv, graupelprv, &
dtp, errmsg, errflg)
!
use machine, only: kind_phys

implicit none

integer, intent(in) :: im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, ntrac, imp_physics, imp_physics_gfdl, imp_physics_thompson
integer, intent(in) :: im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, ntrac
integer, intent(in) :: imp_physics, imp_physics_gfdl, imp_physics_thompson, imp_physics_mg
logical, intent(in) :: cal_pre, lssav, ldiag3d, cplflx, cplchm

real(kind=kind_phys), intent(in) :: dtf, frain, con_g
Expand All @@ -196,8 +201,11 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt
real(kind=kind_phys), dimension(im,levs+1), intent(in) :: prsi, phii
real(kind=kind_phys), dimension(im,levs,ntrac), intent(in) :: gq0

real(kind=kind_phys), dimension(im), intent(inout) :: rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, &
srflag, totprcp, totice, totsnw, totgrp, totprcpb, toticeb, totsnwb, totgrpb, rain_cpl, rainc_cpl, snow_cpl, pwat
real(kind=kind_phys), dimension(im), intent(in ) :: sr
real(kind=kind_phys), dimension(im), intent(inout) :: rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, &
srflag, cnvprcp, totprcp, totice, totsnw, totgrp, cnvprcpb, &
totprcpb, toticeb, totsnwb, totgrpb, rain_cpl, rainc_cpl, &
snow_cpl, pwat
real(kind=kind_phys), dimension(im,levs), intent(inout) :: dt3dt, dq3dt

! Stochastic physics / surface perturbations
Expand All @@ -223,24 +231,22 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt

! DH* TODO: CLEANUP, all of these should be coming in through the argument list
real(kind=kind_phys), parameter :: con_p001= 0.001d0
real(kind=kind_phys), parameter :: con_day = 86400.d0
#ifdef TRANSITION
real(kind=kind_phys), parameter :: con_day = 86400.0d0
real(kind=kind_phys), parameter :: rainmin = 1.0d-13
#else
real(kind=kind_phys), parameter :: rainmin = 1.0e-13
#endif
real(kind=kind_phys), parameter :: p850 = 85000.0
real(kind=kind_phys), parameter :: p850 = 85000.0d0
! *DH

integer :: i, k, ic

real(kind=kind_phys), parameter :: zero = 0.0d0, one = 1.0d0
real(kind=kind_phys) :: crain, csnow, onebg, tem, total_precip
real(kind=kind_phys), dimension(im) :: domr, domzr, domip, doms, t850, work1

! Initialize CCPP error handling variables
errmsg = ''
errflg = 0

onebg = 1.0d0/con_g
onebg = one/con_g

do i = 1, im
rain(i) = rainc(i) + frain * rain1(i) ! time-step convective plus explicit
Expand Down Expand Up @@ -308,6 +314,14 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt
end if
enddo
endif
if (lssav) then
do i=1,im
domr_diag(i) = domr_diag(i) + domr(i) * dtf
domzr_diag(i) = domzr_diag(i) + domzr(i) * dtf
domip_diag(i) = domip_diag(i) + domip(i) * dtf
doms_diag(i) = doms_diag(i) + doms(i) * dtf
enddo
endif

endif

Expand All @@ -316,21 +330,17 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt
! 'totprcpb=', Diag%totprcpb(1),'totprcp=',Diag%totprcp(1), &
! 'rain=',Diag%rain(1)
do i=1,im
cnvprcp (i) = cnvprcp (i) + rainc(i)
totprcp (i) = totprcp (i) + rain(i)
totice (i) = totice (i) + ice(i)
totsnw (i) = totsnw (i) + snow(i)
totgrp (i) = totgrp (i) + graupel(i)

cnvprcpb(i) = cnvprcpb(i) + rainc(i)
totprcpb(i) = totprcpb(i) + rain(i)
toticeb (i) = toticeb (i) + ice(i)
totsnwb (i) = totsnwb (i) + snow(i)
totgrpb (i) = totgrpb (i) + graupel(i)
!
if (cal_pre) then
domr_diag(i) = domr_diag(i) + domr(i) * dtf
domzr_diag(i) = domzr_diag(i) + domzr(i) * dtf
domip_diag(i) = domip_diag(i) + domip(i) * dtf
doms_diag(i) = doms_diag(i) + doms(i) * dtf
endif
enddo

if (ldiag3d) then
Expand All @@ -355,14 +365,16 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt
enddo
enddo

! Conversion factor mm per physics timestep to m per day
tem = dtp * con_p001 / con_day

!> - For GFDL and Thompson MP scheme, determine convective snow by surface temperature;
!! and determine explicit rain/snow by snow/ice/graupel coming out directly from MP
!! and convective rainfall from the cumulus scheme if the surface temperature is below
!! \f$0^oC\f$.
if (imp_physics == imp_physics_gfdl .or. imp_physics == imp_physics_thompson) then
! determine convective rain/snow by surface temperature
! determine large-scale rain/snow by rain/snow coming out directly from MP
tem = dtp * con_p001 / con_day
do i = 1, im
!tprcp(i) = max(0.0, rain(i) )! clu: rain -> tprcp ! DH now lines 245-250
srflag(i) = 0. ! clu: default srflag as 'rain' (i.e. 0)
Expand All @@ -384,22 +396,29 @@ subroutine GFS_MP_generic_post_run(im, ix, levs, kdt, nrcm, ncld, nncl, ntcw, nt
endif
enddo
elseif( .not. cal_pre) then
do i = 1, im
tprcp(i) = max(0.0, rain(i) )! clu: rain -> tprcp
srflag(i) = 0. ! clu: default srflag as 'rain' (i.e. 0)
if (t850(i) <= 273.16) then
srflag(i) = 1. ! clu: set srflag to 'snow' (i.e. 1)
endif
enddo
if (imp_physics == imp_physics_mg) then ! MG microphysics
do i=1,im
if (rain(i)*tem > rainmin) then
srflag(i) = max(zero, min(one, (rain(i)-rainc(i))*sr(i)/rain(i)))
else
srflag(i) = 0.0
endif
enddo
else
do i = 1, im
tprcp(i) = max(0.0, rain(i) )! clu: rain -> tprcp
srflag(i) = 0.0 ! clu: default srflag as 'rain' (i.e. 0)
if (t850(i) <= 273.16) then
srflag(i) = 1.0 ! clu: set srflag to 'snow' (i.e. 1)
endif
enddo
endif
endif

if (cplflx .or. cplchm) then
do i = 1, im
if (t850(i) > 273.16) then
rain_cpl(i) = rain_cpl(i) + rain(i)
else
snow_cpl(i) = snow_cpl(i) + rain(i)
endif
rain_cpl(i) = rain_cpl(i) + rain(i) * (one-srflag(i))
snow_cpl(i) = snow_cpl(i) + rain(i) * srflag(i)
enddo
endif

Expand Down
Loading

0 comments on commit f7915b9

Please sign in to comment.