Skip to content

Commit

Permalink
turn on/off spec_adv option is working in CCPP F-A scheme.
Browse files Browse the repository at this point in the history
  • Loading branch information
mzhangw committed Sep 19, 2019
1 parent 370d49f commit 957ff82
Show file tree
Hide file tree
Showing 5 changed files with 88 additions and 37 deletions.
51 changes: 29 additions & 22 deletions physics/GFS_suite_interstitial.F90
Original file line number Diff line number Diff line change
Expand Up @@ -622,7 +622,7 @@ end subroutine GFS_suite_interstitial_4_finalize
subroutine GFS_suite_interstitial_4_run (im, levs, ltaerosol, lgocart, cplchm, tracers_total, ntrac, ntcw, ntiw, ntclamt, &
ntrw, ntsw, ntrnc, ntsnc, ntgl, ntgnc, ntlnc, ntinc, nn, imp_physics, imp_physics_gfdl, imp_physics_thompson, &
imp_physics_zhao_carr, imp_physics_zhao_carr_pdf,imp_physics_fer_hires, dtf, save_qc, save_qi, con_pi, epsq, &
gq0, clw, cwm, f_ice, f_rain, dqdti, errmsg, errflg)
gq0, clw, cwm, f_ice, f_rain, f_rimef, dqdti,mpirank,mpiroot, errmsg, errflg)

use machine, only: kind_phys

Expand All @@ -646,8 +646,12 @@ subroutine GFS_suite_interstitial_4_run (im, levs, ltaerosol, lgocart, cplchm, t
real(kind=kind_phys), dimension(im,levs), intent(inout) :: cwm
real(kind=kind_phys), dimension(im,levs), intent(inout) :: f_ice
real(kind=kind_phys), dimension(im,levs), intent(inout) :: f_rain
real(kind=kind_phys), dimension(im,levs), intent(inout) :: f_rimef
! dqdti may not be allocated
real(kind=kind_phys), dimension(:,:), intent(inout) :: dqdti
integer, intent(in) :: mpirank
integer, intent(in) :: mpiroot


character(len=*), intent(out) :: errmsg
integer, intent(out) :: errflg
Expand Down Expand Up @@ -720,29 +724,32 @@ subroutine GFS_suite_interstitial_4_run (im, levs, ltaerosol, lgocart, cplchm, t

!MZ
if (imp_physics == imp_physics_fer_hires) then
!MZ: Update CWM,F_ICE,F_RAIN arrays from separate species advection (spec_adv=T)
!MZ: Update CWM,F_ICE,F_RAIN arrays from separate species advection
!(spec_adv=T.or.F)
DO K=1,levs
DO I=1,IM
CWM(I,K)= max(0.0,gq0(i,k,ntcw))+max(0.0,gq0(i,k,ntiw)) &
+max(0.0,gq0(i,k,ntrw))
IF (gq0(I,K,ntiw)>EPSQ) THEN
F_ICE(I,K)=gq0(I,K,ntiw)/CWM(I,K)
ELSE
F_ICE(I,K)=0.0
ENDIF
IF (gq0(I,K,ntrw)>EPSQ) THEN
F_RAIN(I,K)=gq0(I,K,ntrw)/(gq0(I,K,ntcw)+gq0(I,K,ntrw))
ELSE
F_RAIN(I,K)=0.
ENDIF
ENDDO
DO I=1,IM
CWM(I,K)= max(0.0,gq0(i,k,ntcw))+max(0.0,gq0(i,k,ntiw)) &
+max(0.0,gq0(i,k,ntrw))
IF (gq0(I,K,ntiw)>EPSQ) THEN
F_ICE(I,K)=MAX(0.0,MIN(1.,gq0(I,K,ntiw)/CWM(I,K)))
ELSE
F_ICE(I,K)=0.0
ENDIF
IF (gq0(I,K,ntrw)>EPSQ) THEN
F_RAIN(I,K)=gq0(I,K,ntrw)/(gq0(I,K,ntcw)+gq0(I,K,ntrw))
ELSE
F_RAIN(I,K)=0.
ENDIF
ENDDO
ENDDO
!if(mpirank == mpiroot) write (0,*)'interstitial_4: cwm =', &
! maxval(cwm),minval(cwm)
!if(mpirank == mpiroot) write (0,*)'interstitial_4: f_ice =', &
! maxval(f_ice),minval(f_ice)
!if(mpirank == mpiroot) write (0,*)'interstitial_4: f_rain =', &
! maxval(f_rain),minval(f_rain)
if(mpirank == mpiroot) then
write (0,*)'interstitial_4: cwm =', &
maxval(cwm),minval(cwm)
write (0,*)'interstitial_4: f_ice =', &
maxval(f_ice),minval(f_ice)
write (0,*)'interstitial_4: f_rain =', &
maxval(f_rain),minval(f_rain)
end if

endif

Expand Down
25 changes: 25 additions & 0 deletions physics/GFS_suite_interstitial.meta
Original file line number Diff line number Diff line change
Expand Up @@ -1604,6 +1604,15 @@
kind = kind_phys
intent = inout
optional = F
[f_rimef]
standard_name = rime_factor
long_name = rime factor
units = frac
dimensions = (horizontal_dimension,vertical_dimension)
type = real
kind = kind_phys
intent = inout
optional = F
[dtf]
standard_name = time_step_for_dynamics
long_name = dynamics timestep
Expand Down Expand Up @@ -1676,6 +1685,22 @@
kind = kind_phys
intent = inout
optional = F
[mpirank]
standard_name = mpi_rank
long_name = current MPI-rank
units = index
dimensions = ()
type = integer
intent = in
optional = F
[mpiroot]
standard_name = mpi_root
long_name = master MPI-rank
units = index
dimensions = ()
type = integer
intent = in
optional = F
[errmsg]
standard_name = ccpp_error_message
long_name = error message for error handling in CCPP
Expand Down
8 changes: 6 additions & 2 deletions physics/module_MP_FER_HIRES.F90
Original file line number Diff line number Diff line change
Expand Up @@ -240,7 +240,7 @@ SUBROUTINE FER_HIRES (DT,RHgrd, &
& q,qt, &
& LOWLYR,SR, &
& F_ICE_PHY,F_RAIN_PHY,F_RIMEF_PHY, &
& QC,QR,QS, &
& QC,QR,QS, &
& RAINNC,RAINNCV, &
& threads, &
& ims,ime, jms,jme, lm, &
Expand Down Expand Up @@ -611,7 +611,6 @@ SUBROUTINE FER_HIRES (DT,RHgrd, &
QC(I,J,K)=QC(I,J,K)-QR(I,J,K)
ENDIF
ENDIF
!
ENDDO !- i
ENDDO !- k
ENDDO !- j
Expand Down Expand Up @@ -2216,6 +2215,7 @@ SUBROUTINE EGCP01COLUMN_hr ( ARAIN, ASNOW, DTPH, RHC_col, &
!--------- Produces accurate calculation of cloud condensation ---------
!#######################################################################
!
!>\ingroup hafs_famp
REAL FUNCTION CONDENSE (PP,QW,TK,WV,RHgrd,I,J,L)
!
!---------------------------------------------------------------------------------
Expand Down Expand Up @@ -2288,6 +2288,7 @@ END FUNCTION CONDENSE
!---------------- Calculate ice deposition at T<T_ICE ------------------
!#######################################################################
!
!>\ingroup hafs_famp
REAL FUNCTION DEPOSIT (PP,Tdum,WVdum,RHgrd,I,J,L) !-- Debug 20120111
!
!--- Also uses the Asai (1965) algorithm, but uses a different target
Expand Down Expand Up @@ -2448,6 +2449,7 @@ END SUBROUTINE EGCP01COLUMN_hr
! SH 0211/2002

!-----------------------------------------------------------------------
!>\ingroup hafs_famp
SUBROUTINE FERRIER_INIT_hr (GSMDT,MPI_COMM_COMP,MYPE,mpiroot,THREADS)
!-----------------------------------------------------------------------
!-------------------------------------------------------------------------------
Expand Down Expand Up @@ -2741,6 +2743,7 @@ SUBROUTINE FERRIER_INIT_hr (GSMDT,MPI_COMM_COMP,MYPE,mpiroot,THREADS)
!-----------------------------------------------------------------------
END SUBROUTINE FERRIER_INIT_hr
!
!>\ingroup hafs_famp
SUBROUTINE MY_GROWTH_RATES_NMM_hr (DTPH)
!
!--- Below are tabulated values for the predicted mass of ice crystals
Expand Down Expand Up @@ -2786,6 +2789,7 @@ END SUBROUTINE MY_GROWTH_RATES_NMM_hr
!--------- Old GFS saturation vapor pressure lookup tables -----------
!-----------------------------------------------------------------------
!
!>\ingroup hafs_famp
SUBROUTINE GPVS_hr
! ******************************************************************
!$$$ SUBPROGRAM DOCUMENTATION BLOCK
Expand Down
18 changes: 14 additions & 4 deletions physics/module_mp_fer_hires_pre.F90
Original file line number Diff line number Diff line change
Expand Up @@ -31,21 +31,23 @@ end subroutine mp_fer_hires_pre_finalize
!! | cwm | total_cloud_condensate_mixing_ratio_updated_by_physics| total cloud condensate mixing ratio (except water vapor) updated by physics | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | f_ice | fraction_of_ice_water_cloud | fraction of ice water cloud | frac | 2 | real | kind_phys | inout | F |
!! | f_rain | fraction_of_rain_water_cloud | fraction of rain water cloud | frac | 2 | real | kind_phys | inout | F |
!! | f_rimef | rime_factor | rime factor | frac | 2 | real | kind_phys | inout | F |
!! | epsq | minimum_value_of_specific_humidity | floor value for specific humidity | kg kg-1 | 0 | real | kind_phys | in | F |
!! | t | air_temperature | model layer mean temperature | K | 2 | real | kind_phys | inout | F |
!! | qc | cloud_condensed_water_mixing_ratio_updated_by_physics | moist (dry+vapor, no condensates) mixing ratio of cloud condensed water updated by physics | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | qr | rain_water_mixing_ratio_updated_by_physics | moist (dry+vapor, no condensates) mixing ratio of rain water updated by physics | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | qi | ice_water_mixing_ratio_updated_by_physics | moist (dry+vapor, no condensates) mixing ratio of ice water updated by physics | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | qg | mass_weighted_rime_factor_updated_by_physics | mass weighted rime factor updaed by physics | kg kg-1 | 2 | real | kind_phys | inout | F |
!! | spec_adv | flag_for_individual_cloud_species_advected | flag for individual cloud species advected | flag | 0 | logical | | in | F |
!! | kdt | index_of_time_step | current forecast interation | index | 0 | integer | | in | F |
!! | lm | vertical_dimension | number of vertical levels | count | 0 | integer | | in | F |
!! | ime | horizontal_dimension | horizontal dimension | count | 0 | integer | | in | F |
!! | errmsg | ccpp_error_message | error message for error handling in CCPP | none | 0 | character | len=* | out | F |
!! | errflg | ccpp_error_flag | error flag for error handling in CCPP | flag | 0 | integer | | out | F |
!!
subroutine mp_fer_hires_pre_run (CWM,F_ICE,F_RAIN &
subroutine mp_fer_hires_pre_run (CWM,F_ICE,F_RAIN,F_RIMEF &
,EPSQ &
,T,QC,QR,QI &
,T,QC,QR,QI,QG &
,SPEC_ADV,kdt &
,LM,IME,errmsg,errflg )

Expand All @@ -65,8 +67,9 @@ subroutine mp_fer_hires_pre_run (CWM,F_ICE,F_RAIN &
REAL(kind=kind_phys),DIMENSION(1:IME,1:LM),INTENT(INOUT) :: CWM &
,F_ICE &
,F_RAIN &
,F_RIMEF &
,T,QC,QR &
,QI !QG QS
,QI,QG ! QS
!
!--------------------
!-- Local Variables
Expand Down Expand Up @@ -101,15 +104,18 @@ subroutine mp_fer_hires_pre_run (CWM,F_ICE,F_RAIN &
! or at the initial time step
DO K=1,LM
DO I=1,IME
CWM(I,K)=QC(I,K)+QR(I,K)+QI(I,K)
IF (CWM(I,K)>EPSQ) THEN
LIQW=(1.-F_ice(I,K))*CWM(I,K)
QC(I,K)=(1.-F_rain(I,K))*LIQW
QR(I,K)=F_rain(I,K)*LIQW
QI(I,K)=F_ice(I,K)*CWM(I,K)
QG(I,K)=F_RIMEF(I,K)*QI(I,K)
ELSE
QC(I,K)=0.
QR(I,K)=0.
QI(I,K)=0.
QG(I,K)=0.
ENDIF
ENDDO
ENDDO
Expand All @@ -120,9 +126,13 @@ subroutine mp_fer_hires_pre_run (CWM,F_ICE,F_RAIN &
DO I=1,IME
CWM(I,K)=QC(I,K)+QR(I,K)+QI(I,K)
IF (QI(I,K)>EPSQ) THEN
F_ICE(I,K)=QI(I,K)/CWM(I,K)
!F_ICE(I,K)=QI(I,K)/CWM(I,K)
F_ICE(I,K) = MAX(0., MIN(1., QI(I,K)/CWM(I,K)))
!MZ: f_rimef= qrimef/qi
F_RIMEF(I,K) = QG(I,K)/QI(I,K)
ELSE
F_ICE(I,K)=0.0
F_RIMEF(I,K) = 1.0
ENDIF
IF (QR(I,K)>EPSQ) THEN
F_RAIN(I,K)=QR(I,K)/(QC(I,K)+QR(I,K))
Expand Down
Loading

0 comments on commit 957ff82

Please sign in to comment.