-
Notifications
You must be signed in to change notification settings - Fork 153
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Radar-derived microphysics temperature tendencies similar to operational HRRR #823
Changes from all commits
f1a00a6
914482a
c213000
e687c04
037e564
d247e33
7072fbe
64d8f21
2af3ca3
db934c8
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -13,12 +13,12 @@ end subroutine GFS_MP_generic_pre_init | |
!! \htmlinclude GFS_MP_generic_pre_run.html | ||
!! | ||
subroutine GFS_MP_generic_pre_run(im, levs, ldiag3d, qdiag3d, do_aw, ntcw, nncl, & | ||
ntrac, gt0, gq0, save_t, save_q, errmsg, errflg) | ||
ntrac, gt0, gq0, save_t, save_q, num_dfi_radar, errmsg, errflg) | ||
! | ||
use machine, only: kind_phys | ||
|
||
implicit none | ||
integer, intent(in) :: im, levs, ntcw, nncl, ntrac | ||
integer, intent(in) :: im, levs, ntcw, nncl, ntrac, num_dfi_radar | ||
logical, intent(in) :: ldiag3d, qdiag3d, do_aw | ||
real(kind=kind_phys), dimension(:,:), intent(in) :: gt0 | ||
real(kind=kind_phys), dimension(:,:,:), intent(in) :: gq0 | ||
|
@@ -35,12 +35,14 @@ subroutine GFS_MP_generic_pre_run(im, levs, ldiag3d, qdiag3d, do_aw, ntcw, nncl, | |
errmsg = '' | ||
errflg = 0 | ||
|
||
if (ldiag3d .or. do_aw) then | ||
if (ldiag3d .or. do_aw .or. num_dfi_radar>0) then | ||
do k=1,levs | ||
do i=1,im | ||
save_t(i,k) = gt0(i,k) | ||
enddo | ||
enddo | ||
endif | ||
if (ldiag3d .or. do_aw) then | ||
if(qdiag3d) then | ||
do n=1,ntrac | ||
do k=1,levs | ||
|
@@ -91,28 +93,36 @@ subroutine GFS_MP_generic_post_run( | |
graupel0, del, rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, srflag, sr, cnvprcp, totprcp, totice, & | ||
totsnw, totgrp, cnvprcpb, totprcpb, toticeb, totsnwb, totgrpb, rain_cpl, rainc_cpl, snow_cpl, pwat, & | ||
drain_cpl, dsnow_cpl, lsm, lsm_ruc, lsm_noahmp, raincprv, rainncprv, iceprv, snowprv, & | ||
graupelprv, draincprv, drainncprv, diceprv, dsnowprv, dgraupelprv, dtp, & | ||
dtend, dtidx, index_of_temperature, index_of_process_mp,ldiag3d, qdiag3d, lssav, & | ||
errmsg, errflg) | ||
graupelprv, draincprv, drainncprv, diceprv, dsnowprv, dgraupelprv, dtp, dfi_radar_max_intervals, & | ||
dtend, dtidx, index_of_temperature, index_of_process_mp,ldiag3d, qdiag3d, lssav, num_dfi_radar, fh_dfi_radar, & | ||
index_of_process_dfi_radar, ix_dfi_radar, dfi_radar_tten, radar_tten_limits, fhour, errmsg, errflg) | ||
! | ||
use machine, only: kind_phys | ||
|
||
implicit none | ||
|
||
integer, intent(in) :: im, levs, kdt, nrcm, nncl, ntcw, ntrac | ||
integer, intent(in) :: im, levs, kdt, nrcm, nncl, ntcw, ntrac, num_dfi_radar, index_of_process_dfi_radar | ||
integer, intent(in) :: imp_physics, imp_physics_gfdl, imp_physics_thompson, imp_physics_mg, imp_physics_fer_hires | ||
logical, intent(in) :: cal_pre, lssav, ldiag3d, qdiag3d, cplflx, cplchm | ||
integer, intent(in) :: index_of_temperature,index_of_process_mp | ||
|
||
integer :: dfi_radar_max_intervals | ||
real(kind=kind_phys), intent(in) :: fh_dfi_radar(:), fhour | ||
real(kind=kind_phys), intent(in) :: radar_tten_limits(:) | ||
integer :: ix_dfi_radar(:) | ||
real(kind=kind_phys), dimension(:,:), intent(inout) :: gt0 | ||
|
||
real(kind=kind_phys), intent(in) :: dtf, frain, con_g, rainmin | ||
real(kind=kind_phys), dimension(:), intent(in) :: rain1, xlat, xlon, tsfc | ||
real(kind=kind_phys), dimension(:), intent(inout) :: ice, snow, graupel, rainc | ||
real(kind=kind_phys), dimension(:), intent(in) :: rain0, ice0, snow0, graupel0 | ||
real(kind=kind_phys), dimension(:,:), intent(in) :: rann | ||
real(kind=kind_phys), dimension(:,:), intent(in) :: gt0, prsl, save_t, del | ||
real(kind=kind_phys), dimension(:,:), intent(in) :: prsl, save_t, del | ||
real(kind=kind_phys), dimension(:,:), intent(in) :: prsi, phii | ||
real(kind=kind_phys), dimension(:,:,:), intent(in) :: gq0, save_q | ||
|
||
real(kind=kind_phys), dimension(:,:,:), intent(in) :: dfi_radar_tten | ||
|
||
real(kind=kind_phys), dimension(:), intent(in ) :: sr | ||
real(kind=kind_phys), dimension(:), intent(inout) :: rain, domr_diag, domzr_diag, domip_diag, doms_diag, tprcp, & | ||
srflag, cnvprcp, totprcp, totice, totsnw, totgrp, cnvprcpb, & | ||
|
@@ -150,10 +160,10 @@ subroutine GFS_MP_generic_post_run( | |
real(kind=kind_phys), parameter :: p850 = 85000.0_kind_phys | ||
! *DH | ||
|
||
integer :: i, k, ic, itrac, idtend | ||
integer :: i, k, ic, itrac, idtend, itime, idtend_radar, idtend_mp | ||
|
||
real(kind=kind_phys), parameter :: zero = 0.0_kind_phys, one = 1.0_kind_phys | ||
real(kind=kind_phys) :: crain, csnow, onebg, tem, total_precip, tem1, tem2 | ||
real(kind=kind_phys) :: crain, csnow, onebg, tem, total_precip, tem1, tem2, ttend | ||
real(kind=kind_phys), dimension(im) :: domr, domzr, domip, doms, t850, work1 | ||
|
||
! Initialize CCPP error handling variables | ||
|
@@ -244,6 +254,52 @@ subroutine GFS_MP_generic_post_run( | |
|
||
endif | ||
|
||
do itime=1,num_dfi_radar | ||
if(ix_dfi_radar(itime)<1) cycle | ||
if(fhour<fh_dfi_radar(itime)) cycle | ||
if(fhour>=fh_dfi_radar(itime+1)) cycle | ||
exit | ||
enddo | ||
if_radar: if(itime<=num_dfi_radar) then | ||
radar_k: do k=3,levs-2 ! Avoid model top and bottom in case DA forgets to | ||
radar_i: do i=1,im | ||
ttend = dfi_radar_tten(i,k,itime) | ||
if_active: if (ttend>-19) then | ||
ttend = max(ttend,radar_tten_limits(1)) | ||
ttend = min(ttend,radar_tten_limits(2)) | ||
|
||
! add radar temp tendency | ||
! there is radar coverage | ||
gt0(i,k) = save_t(i,k) + ttend*dtp | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I have reservations about this functionality being added to this interstitial scheme. Other than one interstitial scheme that updates the state variables after the process-split section of physics is completed, I don't think that any others update the state. In general, they are supposed to have "glue code", not necessarily algorithms that actually affect the state variables. Why could this functionality not be added to a new scheme, especially if it is supposed to be used with HRRR-based suites at this point? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. In other words, this is sorta a physics parameterization itself, "helping out" microphysics. Why should it not exist as its own scheme? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. It can't because the calcpreciptype has to use the actual microphysics temperature tendencies, or it'll get meaningless values. Somehow it has to go after that calculation, but before the later code. There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. OK, I see what you mean. It's a bit messy and goes against precedent, but given how the code is organized now, I don't see an alternative if the constraint you mentioned is legit. |
||
end if if_active | ||
end do radar_i | ||
end do radar_k | ||
if(ldiag3d) then | ||
idtend_radar = dtidx(index_of_temperature,index_of_process_dfi_radar) | ||
idtend_mp = dtidx(index_of_temperature,index_of_process_mp) | ||
if(idtend_radar>0 .or. idtend_mp>0) then | ||
if(idtend_mp>0) then | ||
dtend(:,1:2,idtend_mp) = dtend(:,1:2,idtend_mp) + (gt0(:,1:2)-save_t(:,1:2))*frain | ||
endif | ||
do k=3,levs-2 ! Avoid model top and bottom in case DA forgets to | ||
do i=1,im | ||
ttend = dfi_radar_tten(i,k,itime) | ||
if (ttend>-19) then | ||
if(idtend_radar>0) then | ||
dtend(i,k,idtend_radar) = dtend(i,k,idtend_radar) + (gt0(i,k)-save_t(i,k)) * frain | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Aren't you saving microphysics tendency + radar tendency in the radar tendency slice here? Won't it be possible for both microphysics and radar tendencies to be nonzero simultaneously? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. The radar tendencies override the microphysics tendency; that is the purpose of this PR. In each gridpoint, you can have one, or the other, not both. There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. At the point where this interstitial scheme is called, microphysics tendencies have already been applied since they are being used as time-split in every CCPP suite that I'm aware of. Unless I'm mistaken, this is only replacing the microphysics tendency diagnostic, not the one that actually gets applied to the state. Wouldn't one need to insert code in every microphysics scheme to bypass it when you wanted to apply the radar tendencies? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. No. In the RAP and HRRR suites, only the thompson scheme is run between GFS_MP_generic_pre and _post. There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. OK, so the Thompson scheme is called and updates gt0. Then, this is called and replaces the previously calculated gt0 if a radar tendency is present. |
||
endif | ||
else if(idtend_mp>0) then | ||
dtend(i,k,idtend_mp) = dtend(i,k,idtend_mp) + (gt0(i,k)-save_t(i,k)) * frain | ||
endif | ||
enddo | ||
enddo | ||
if(idtend_mp>0) then | ||
dtend(:,levs-1:levs,idtend_mp) = dtend(:,levs-1:levs,idtend_mp) + (gt0(:,levs-1:levs)-save_t(:,levs-1:levs))*frain | ||
endif | ||
endif | ||
endif | ||
endif if_radar | ||
|
||
t850(1:im) = gt0(1:im,1) | ||
|
||
do k = 1, levs-1 | ||
|
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I don't understand the hard-coded test for >-19. What if a user specifies radar_tten_limits(1) < -19?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Such a large value would be meaningless. The fill value in the data assimilation code is -20, and that is what the "if" block tests for.