Merge remote-tracking branch 'origin/release-v4.6.1'
[WRF.git] / phys / module_cam_mp_conv_water.F
bloba0503b18769950305f3dd3ad544cb8ec97d68a4a
1 #define WRF_PORT
2 #define MODAL_AERO
3   module conv_water
5    ! --------------------------------------------------------------------- ! 
6    ! Purpose:                                                              !
7    ! Computes grid-box average liquid (and ice) from stratus and cumulus   !
8    ! Just for the purposes of radiation.                                   !
9    !                                                                       ! 
10    ! Method:                                                               !
11    ! Extract information about deep+shallow liquid and cloud fraction from !
12    ! the physics buffer.                                                   !
13    !                                                                       !
14    ! Author: Rich Neale, August 2006                                       !
15    !         October 2006: Allow averaging of liquid to give a linear      !
16    !                       average in emissivity.                          !
17    !         Andrew Gettelman October 2010  Separate module                !
18    !---------------------------------------------------------------------- !
20   use shr_kind_mod,  only: r8=>shr_kind_r8
21 #ifndef WRF_PORT
22   use ppgrid,        only: pcols, pver, pverp
23 #else
24   use module_cam_support, only: pcols, pver, pverp
25 #endif
26   use physconst,     only: gravit, latvap, latice
27 #ifndef WRF_PORT
28   use abortutils,    only: endrun
30   use perf_mod
31   use cam_logfile,   only: iulog
32 #else
33   use module_cam_support, only: endrun, iulog
34 #endif
36   implicit none
37   private
38   save
40   public :: conv_water_register, conv_water_4rad, conv_water_init
42 ! pbuf indices
44   integer :: icwmrsh_idx, icwmrdp_idx, fice_idx, sh_frac_idx, dp_frac_idx, concldql_idx, &
45              ast_idx, alst_idx, aist_idx, qlst_idx, qist_idx, sh_cldliq1_idx, sh_cldice1_idx
47   contains
49   !============================================================================ !
51   subroutine conv_water_register
53   !---------------------------------------------------------------------- !
54   !                                                                       !
55   ! Register the fields in the physics buffer.                            !
56   !                                                                       !
57   !---------------------------------------------------------------------- !
58 #ifndef WRF_PORT
59     use constituents, only: cnst_add, pcnst
60     use physconst,    only: mwdry, cpair
61     use phys_buffer,  only: pbuf_times, pbuf_add
63   !-----------------------------------------------------------------------
65     ! these calls were already done in convect_shallow...so here I add the same fields to the physics buffer with a "1" at the end
66     call pbuf_add('SH_CLDLIQ1', 'physpkg', 1, pver,  1, sh_cldliq1_idx)  ! shallow gbm cloud liquid water (kg/kg)
67     call pbuf_add('SH_CLDICE1', 'physpkg', 1, pver,  1, sh_cldice1_idx)  ! shallow gbm cloud ice water (kg/kg)
68 #endif
70   end subroutine conv_water_register
73   !============================================================================ !
74   !                                                                             !
75   !============================================================================ !
77    subroutine conv_water_init()
78    ! --------------------------------------------------------------------- ! 
79    ! Purpose:                                                              !
80    !   Initializes the pbuf indices required by conv_water
81    ! --------------------------------------------------------------------- ! 
82 #ifndef WRF_PORT
83    use phys_buffer,     only: pbuf_size_max, pbuf_fld, pbuf_old_tim_idx, pbuf_get_fld_idx 
84 #endif
85    implicit none
86 #ifndef WRF_PORT
87    icwmrsh_idx  = pbuf_get_fld_idx('ICWMRSH')
88    icwmrdp_idx  = pbuf_get_fld_idx('ICWMRDP')
89    fice_idx     = pbuf_get_fld_idx('FICE')
90    sh_frac_idx  = pbuf_get_fld_idx('SH_FRAC')
91    dp_frac_idx  = pbuf_get_fld_idx('DP_FRAC')
92    concldql_idx = pbuf_get_fld_idx('CONCLDQL')
93    ast_idx      = pbuf_get_fld_idx('AST')
94    alst_idx     = pbuf_get_fld_idx('ALST')
95    aist_idx     = pbuf_get_fld_idx('AIST')
96    qlst_idx     = pbuf_get_fld_idx('QLST')
97    qist_idx     = pbuf_get_fld_idx('QIST')
98 #endif
99    end subroutine conv_water_init
100 #ifndef WRF_PORT
101    subroutine conv_water_4rad( lchnk, ncol, pbuf, conv_water_mode, &
102                                rei, pdel, ls_liq, ls_ice, totg_liq, totg_ice )
103 #else
104      !Replace pbuf with actual variables
105    subroutine conv_water_4rad( lchnk, ncol, ast, sh_icwmr, dp_icwmr, &
106         fice, sh_frac, dp_frac, conv_water_mode, rei, pdel, ls_liq,  &
107         ls_ice, totg_liq, totg_ice )
108 #endif
110    ! --------------------------------------------------------------------- ! 
111    ! Purpose:                                                              !
112    ! Computes grid-box average liquid (and ice) from stratus and cumulus   !
113    ! Just for the purposes of radiation.                                   !
114    !                                                                       ! 
115    ! Method:                                                               !
116    ! Extract information about deep+shallow liquid and cloud fraction from !
117    ! the physics buffer.                                                   !
118    !                                                                       !
119    ! Author: Rich Neale, August 2006                                       !
120    !         October 2006: Allow averaging of liquid to give a linear      !
121    !                       average in emissivity.                          !
122    !                                                                       !
123    !---------------------------------------------------------------------- !
124 #ifndef WRF_PORT
125    use phys_buffer,     only: pbuf_size_max, pbuf_fld, pbuf_old_tim_idx, pbuf_get_fld_idx 
126    use cam_history,     only: outfld
127    use phys_control,    only: phys_getopts
128    use phys_debug_util, only: phys_debug_col
129 #else
130    use module_cam_support, only: outfld
131 #endif
132    
133    implicit none
135    ! ---------------------- !
136    ! Input-Output Arguments !
137    ! ---------------------- !
138 #ifndef WRF_PORT
139    type(pbuf_fld), intent(inout), dimension(pbuf_size_max) :: pbuf
140 #endif
141    integer,  intent(in) :: lchnk
142    integer,  intent(in) :: ncol
143    integer,  intent(in) :: conv_water_mode
144    real(r8), intent(in) :: rei(pcols,pver)        ! Ice effective drop size (microns)
145    real(r8), intent(in) :: pdel(pcols,pver)       ! Moist pressure difference across layer
146    real(r8), intent(in) :: ls_liq(pcols,pver)     ! Large-scale contributions to GBA cloud liq      
147    real(r8), intent(in) :: ls_ice(pcols,pver)     ! Large-scale contributions to GBA cloud ice 
148 #ifdef WRF_PORT
149    real(r8), intent(in) :: ast(pcols,pver)
150    real(r8), intent(in) :: sh_icwmr(pcols,pver)
151    real(r8), intent(in) :: dp_icwmr(pcols,pver)
152    real(r8), intent(in) :: fice(pcols,pver)
153    real(r8), intent(in) :: sh_frac(pcols,pver)
154    real(r8), intent(in) :: dp_frac(pcols,pver)
155 #endif
156    real(r8), intent(out):: totg_ice(pcols,pver)   ! Total GBA in-cloud ice
157    real(r8), intent(out):: totg_liq(pcols,pver)   ! Total GBA in-cloud liquid
159    ! --------------- !
160    ! Local Workspace !
161    ! --------------- !
163    ! Physics buffer fields
164 #ifndef WRF_PORT
165    real(r8), pointer, dimension(:,:) ::  ast      ! Physical liquid+ice stratus cloud fraction
166    real(r8), pointer, dimension(:,:) ::  cu_frac  ! Final convective cloud fraction
167    real(r8), pointer, dimension(:,:) ::  sh_frac  ! Shallow convective cloud fraction
168    real(r8), pointer, dimension(:,:) ::  dp_frac  ! Deep convective cloud fraction
170    real(r8), pointer, dimension(:,:) ::  alst     ! Physical liquid stratus cloud fraction
171    real(r8), pointer, dimension(:,:) ::  aist     ! Physical ice    stratus cloud fraction
172    real(r8), pointer, dimension(:,:) ::  qlst     ! Physical in-stratus LWC [kg/kg]
173    real(r8), pointer, dimension(:,:) ::  qist     ! Physical in-stratus IWC [kg/kg]
175    real(r8), pointer, dimension(:,:) ::  dp_icwmr ! Deep conv. cloud water
176    real(r8), pointer, dimension(:,:) ::  sh_icwmr ! Shallow conv. cloud water
177    real(r8), pointer, dimension(:,:) ::  fice     ! Ice partitioning ratio
178    real(r8), pointer, dimension(:,:) ::  sh_cldliq ! shallow convection gbx liq cld mixing ratio for COSP
179    real(r8), pointer, dimension(:,:) ::  sh_cldice ! shallow convection gbx ice cld mixing ratio for COSP
180 #else
181    real(r8), dimension(pcols,pver) ::  sh_cldliq ! shallow convection gbx liq cld mixing ratio for COSP
182    real(r8), dimension(pcols,pver) ::  sh_cldice ! shallow convection gbx ice cld mixing ratio for COSP
183 #endif
185    ! Local Variables
187    real(r8) :: conv_ice(pcols,pver)               ! Convective contributions to IC cloud ice
188    real(r8) :: conv_liq(pcols,pver)               ! Convective contributions to IC cloud liquid
189    real(r8) :: tot_ice(pcols,pver)                ! Total IC ice
190    real(r8) :: tot_liq(pcols,pver)                ! Total IC liquid
192    integer  :: i,k,itim                           ! Lon, lev indices buff stuff.
193    real(r8) :: cu_icwmr                           ! Convective  water for this grid-box.   
194    real(r8) :: ls_icwmr                           ! Large-scale water for this grid-box. 
195    real(r8) :: tot_icwmr                          ! Large-scale water for this grid-box.  
196    real(r8) :: ls_frac                            ! Large-scale cloud frac for this grid-box. 
197    real(r8) :: tot0_frac, cu0_frac, dp0_frac, sh0_frac 
198    real(r8) :: kabs, kabsi, kabsl, alpha, dp0, sh0, ic_limit, frac_limit  
199    real(r8) :: wrk1         
201    ! --------- !
202    ! Parameter !
203    ! --------- !
205    parameter( kabsl = 0.090361_r8, frac_limit = 0.01_r8, ic_limit = 1.e-12_r8 )
207  ! Get microphysics option
209    character(len=16) :: microp_scheme 
210 #ifndef WRF_PORT   
211    call phys_getopts( microp_scheme_out = microp_scheme )
212 #else
213    microp_scheme = 'MG'
214 #endif
216  ! Get convective in-cloud water and ice/water temperature partitioning.
217 #ifndef WRF_PORT   
218    sh_icwmr => pbuf(icwmrsh_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
219    dp_icwmr => pbuf(icwmrdp_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
220    fice => pbuf(fice_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
222  ! Get convective in-cloud fraction    
224    sh_frac => pbuf(sh_frac_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
225    dp_frac => pbuf(dp_frac_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
226    cu_frac => pbuf(concldql_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)  
228    itim = pbuf_old_tim_idx()
229    ast => pbuf(ast_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,itim) 
231    itim = pbuf_old_tim_idx()
232    alst => pbuf(alst_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,itim) 
233    itim = pbuf_old_tim_idx()
234    aist => pbuf(aist_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,itim) 
235    itim = pbuf_old_tim_idx()
236    qlst => pbuf(qlst_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,itim) 
237    itim = pbuf_old_tim_idx()
238    qist => pbuf(qist_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,itim) 
239 #endif
241    ! --------------------------------------------------------------- !
242    ! Loop through grid-boxes and determine:                          !
243    ! 1. Effective mean in-cloud convective ice/liquid (deep+shallow) !
244    ! 2. Effective mean in-cloud total ice/liquid (ls+convective)     !
245    ! --------------------------------------------------------------- !
247    do k = 1, pver
248    do i = 1, ncol
250       if( sh_frac(i,k) <= frac_limit .or. sh_icwmr(i,k) <= ic_limit ) then
251           sh0_frac = 0._r8
252       else
253           sh0_frac = sh_frac(i,k)
254       endif
255       if( dp_frac(i,k) <= frac_limit .or. dp_icwmr(i,k) <= ic_limit ) then
256           dp0_frac = 0._r8
257       else
258           dp0_frac = dp_frac(i,k)
259       endif
260       cu0_frac = sh0_frac + dp0_frac
262     ! For the moment calculate the emissivity based upon the ls clouds ice fraction
264       wrk1 = min(1._r8,max(0._r8, ls_ice(i,k)/(ls_ice(i,k)+ls_liq(i,k)+1.e-36_r8)))
266       if( ( cu0_frac < frac_limit ) .or. ( ( sh_icwmr(i,k) + dp_icwmr(i,k) ) < ic_limit ) ) then
268             cu0_frac = 0._r8
269             cu_icwmr = 0._r8
270          
271             ls_frac = ast(i,k)
272             if( ls_frac < frac_limit ) then
273                 ls_frac  = 0._r8
274                 ls_icwmr = 0._r8
275             else
276                 ls_icwmr = ( ls_liq(i,k) + ls_ice(i,k) )/max(frac_limit,ls_frac) ! Convert to IC value.
277             end if
279             tot0_frac = ls_frac
280             tot_icwmr = ls_icwmr
281            
282       else
284           ! Select radiation constants (effective radii) for emissivity averaging.
285             
286             if( microp_scheme .eq. 'MG' ) then
287                 kabsi = 0.005_r8 + 1._r8/min(max(13._r8,rei(i,k)),130._r8)
288             elseif( microp_scheme .eq. 'RK' ) then
289                 kabsi = 0.005_r8 + 1._r8/rei(i,k)
290             endif
291             kabs  = kabsl * ( 1._r8 - wrk1 ) + kabsi * wrk1
292             alpha = -1.66_r8*kabs*pdel(i,k)/gravit*1000.0_r8
294           ! Selecting cumulus in-cloud water.            
296             select case (conv_water_mode) ! Type of average
297             case (1) ! Area weighted arithmetic average
298                cu_icwmr = ( sh0_frac * sh_icwmr(i,k) + dp0_frac*dp_icwmr(i,k))/max(frac_limit,cu0_frac)
299             case (2)
300                sh0 = exp(alpha*sh_icwmr(i,k))
301                dp0 = exp(alpha*dp_icwmr(i,k))               
302                cu_icwmr = log((sh0_frac*sh0+dp0_frac*dp0)/max(frac_limit,cu0_frac))
303                cu_icwmr = cu_icwmr/alpha
304             case default ! Area weighted 'arithmetic in emissivity' average.
305 !               call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
306             end select
308           ! Selecting total in-cloud water. 
309           ! Attribute large-scale/convective area fraction differently from default.
311             ls_frac   = ast(i,k) 
312             ls_icwmr  = (ls_liq(i,k) + ls_ice(i,k))/max(frac_limit,ls_frac) ! Convert to IC value.
313             tot0_frac = (ls_frac + cu0_frac) 
315             select case (conv_water_mode) ! Type of average
316             case (1) ! Area weighted 'arithmetic in emissivity' average
317                tot_icwmr = (ls_frac*ls_icwmr + cu0_frac*cu_icwmr)/max(frac_limit,tot0_frac)
318             case (2)
319                tot_icwmr = log((ls_frac*exp(alpha*ls_icwmr)+cu0_frac*exp(alpha*cu_icwmr))/max(frac_limit,tot0_frac))
320                tot_icwmr = tot_icwmr/alpha
321             case default ! Area weighted 'arithmetic in emissivity' average.
322 !               call endrun ('CONV_WATER_4_RAD: Unknown option for conv_water_in_rad - exiting')
323             end select
325       end if
327     ! Repartition convective cloud water into liquid and ice phase.
328     ! Currently, this partition is made using the ice fraction of stratus condensate.
329     ! In future, we should use ice fraction explicitly computed from the convection scheme.
331       conv_ice(i,k) = cu_icwmr * wrk1
332       conv_liq(i,k) = cu_icwmr * (1._r8-wrk1)
334       tot_ice(i,k)  = tot_icwmr * wrk1
335       tot_liq(i,k)  = tot_icwmr * (1._r8-wrk1)
337       totg_ice(i,k) = tot0_frac * tot_icwmr * wrk1
338       totg_liq(i,k) = tot0_frac * tot_icwmr * (1._r8-wrk1)
340    end do
341    end do
343 !add pbuff calls for COSP
344 #ifndef WRF_PORT   
345    sh_cldliq  => pbuf(sh_cldliq1_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
346    sh_cldice  => pbuf(sh_cldice1_idx)%fld_ptr(1,1:pcols,1:pver,lchnk,1)
347 #endif
348    sh_cldliq(:ncol,:pver)= sh_icwmr(:ncol,:pver)*(1-fice(:ncol,:pver))*sh_frac(:ncol,:pver)
349    sh_cldice(:ncol,:pver)=sh_icwmr(:ncol,:pver)*fice(:ncol,:pver)*sh_frac(:ncol,:pver)
351   ! Output convective IC WMRs
352    
353    call outfld( 'ICLMRCU ', conv_liq  , pcols, lchnk )
354    call outfld( 'ICIMRCU ', conv_ice  , pcols, lchnk )
355    call outfld( 'ICWMRSH ', sh_icwmr  , pcols, lchnk )
356    call outfld( 'ICWMRDP ', dp_icwmr  , pcols, lchnk ) 
357    call outfld( 'ICLMRTOT', tot_liq   , pcols, lchnk )
358    call outfld( 'ICIMRTOT', tot_ice   , pcols, lchnk )
359    call outfld( 'SH_CLD  ', sh_frac   , pcols, lchnk )
360    call outfld( 'DP_CLD  ', dp_frac   , pcols, lchnk )
362   end subroutine conv_water_4rad
364 end module conv_water