Line data Source code
1 : module pbl_utils
2 : !-----------------------------------------------------------------------!
3 : ! Module to hold PBL-related subprograms that may be used with multiple !
4 : ! different vertical diffusion schemes. !
5 : ! !
6 : ! Public subroutines: !
7 : !
8 : ! calc_obklen !
9 : ! !
10 : !------------------ History --------------------------------------------!
11 : ! Created: Apr. 2012, by S. Santos !
12 : !-----------------------------------------------------------------------!
13 :
14 : use shr_kind_mod, only: r8 => shr_kind_r8
15 :
16 : implicit none
17 : private
18 :
19 : ! Public Procedures
20 : !----------------------------------------------------------------------!
21 : ! Excepting the initialization procedure, these are elemental
22 : ! procedures, so they can accept scalars or any dimension of array as
23 : ! arguments, as long as all arguments have the same number of
24 : ! elements.
25 : public compute_radf
26 :
27 : contains
28 :
29 0 : subroutine compute_radf( choice_radf, i, pcols, pver, ncvmax, ncvfin, ktop, qmin, &
30 0 : ql, pi, qrlw, g, cldeff, zi, chs, lwp_CL, opt_depth_CL, &
31 0 : radinvfrac_CL, radf_CL )
32 : ! -------------------------------------------------------------------------- !
33 : ! Purpose: !
34 : ! Calculate cloud-top radiative cooling contribution to buoyancy production. !
35 : ! Here, 'radf' [m2/s3] is additional buoyancy flux at the CL top interface !
36 : ! associated with cloud-top LW cooling being mainly concentrated near the CL !
37 : ! top interface ( just below CL top interface ). Contribution of SW heating !
38 : ! within the cloud is not included in this radiative buoyancy production !
39 : ! since SW heating is more broadly distributed throughout the CL top layer. !
40 : ! -------------------------------------------------------------------------- !
41 :
42 : !-----------------!
43 : ! Input variables !
44 : !-----------------!
45 : character(len=6), intent(in) :: choice_radf ! Method for calculating radf
46 : integer, intent(in) :: i ! Index of current column
47 : integer, intent(in) :: pcols ! Number of atmospheric columns
48 : integer, intent(in) :: pver ! Number of atmospheric layers
49 : integer, intent(in) :: ncvmax ! Max numbers of CLs (perhaps equal to pver)
50 : integer, intent(in) :: ncvfin(pcols) ! Total number of CL in column
51 : integer, intent(in) :: ktop(pcols, ncvmax) ! ktop for current column
52 : real(r8), intent(in) :: qmin ! Minimum grid-mean LWC counted as clouds [kg/kg]
53 : real(r8), intent(in) :: ql(pcols, pver) ! Liquid water specific humidity [ kg/kg ]
54 : real(r8), intent(in) :: pi(pcols, pver+1) ! Interface pressures [ Pa ]
55 : real(r8), intent(in) :: qrlw(pcols, pver) ! Input grid-mean LW heating rate : [ K/s ] * cpair * dp = [ W/kg*Pa ]
56 : real(r8), intent(in) :: g ! Gravitational acceleration
57 : real(r8), intent(in) :: cldeff(pcols,pver) ! Effective Cloud Fraction [fraction]
58 : real(r8), intent(in) :: zi(pcols, pver+1) ! Interface heights [ m ]
59 : real(r8), intent(in) :: chs(pcols, pver+1) ! Buoyancy coeffi. saturated sl (heat) coef. at all interfaces.
60 :
61 : !------------------!
62 : ! Output variables !
63 : !------------------!
64 : real(r8), intent(out) :: lwp_CL(ncvmax) ! LWP in the CL top layer [ kg/m2 ]
65 : real(r8), intent(out) :: opt_depth_CL(ncvmax) ! Optical depth of the CL top layer
66 : real(r8), intent(out) :: radinvfrac_CL(ncvmax) ! Fraction of LW radiative cooling confined in the top portion of CL
67 : real(r8), intent(out) :: radf_CL(ncvmax) ! Buoyancy production at the CL top due to radiative cooling [ m2/s3 ]
68 :
69 : !-----------------!
70 : ! Local variables !
71 : !-----------------!
72 : integer :: kt, ncv
73 : real(r8) :: lwp, opt_depth, radinvfrac, radf
74 :
75 :
76 : !-----------------!
77 : ! Begin main code !
78 : !-----------------!
79 0 : lwp_CL = 0._r8
80 0 : opt_depth_CL = 0._r8
81 0 : radinvfrac_CL = 0._r8
82 0 : radf_CL = 0._r8
83 :
84 : ! ---------------------------------------- !
85 : ! Perform do loop for individual CL regime !
86 : ! ---------------------------------------- !
87 0 : do ncv = 1, ncvfin(i)
88 0 : kt = ktop(i,ncv)
89 : !-----------------------------------------------------!
90 : ! Compute radf for each CL regime and for each column !
91 : !-----------------------------------------------------!
92 0 : if( choice_radf .eq. 'orig' ) then
93 0 : if( ql(i,kt) .gt. qmin .and. ql(i,kt-1) .lt. qmin ) then
94 0 : lwp = ql(i,kt) * ( pi(i,kt+1) - pi(i,kt) ) / g
95 0 : opt_depth = 156._r8 * lwp ! Estimated LW optical depth in the CL top layer
96 : ! Approximate LW cooling fraction concentrated at the inversion by using
97 : ! polynomial approx to exact formula 1-2/opt_depth+2/(exp(opt_depth)-1))
98 :
99 0 : radinvfrac = opt_depth * ( 4._r8 + opt_depth ) / ( 6._r8 * ( 4._r8 + opt_depth ) + opt_depth**2 )
100 0 : radf = qrlw(i,kt) / ( pi(i,kt) - pi(i,kt+1) ) ! Cp*radiative cooling = [ W/kg ]
101 0 : radf = max( radinvfrac * radf * ( zi(i,kt) - zi(i,kt+1) ), 0._r8 ) * chs(i,kt)
102 : ! We can disable cloud LW cooling contribution to turbulence by uncommenting:
103 : ! radf = 0._r8
104 : end if
105 :
106 0 : elseif( choice_radf .eq. 'ramp' ) then
107 :
108 0 : lwp = ql(i,kt) * ( pi(i,kt+1) - pi(i,kt) ) / g
109 0 : opt_depth = 156._r8 * lwp ! Estimated LW optical depth in the CL top layer
110 0 : radinvfrac = opt_depth * ( 4._r8 + opt_depth ) / ( 6._r8 * ( 4._r8 + opt_depth ) + opt_depth**2 )
111 0 : radinvfrac = max(cldeff(i,kt)-cldeff(i,kt-1),0._r8) * radinvfrac
112 0 : radf = qrlw(i,kt) / ( pi(i,kt) - pi(i,kt+1) ) ! Cp*radiative cooling [W/kg]
113 0 : radf = max( radinvfrac * radf * ( zi(i,kt) - zi(i,kt+1) ), 0._r8 ) * chs(i,kt)
114 :
115 0 : elseif( choice_radf .eq. 'maxi' ) then
116 :
117 : ! Radiative flux divergence both in 'kt' and 'kt-1' layers are included
118 : ! 1. From 'kt' layer
119 0 : lwp = ql(i,kt) * ( pi(i,kt+1) - pi(i,kt) ) / g
120 0 : opt_depth = 156._r8 * lwp ! Estimated LW optical depth in the CL top layer
121 0 : radinvfrac = opt_depth * ( 4._r8 + opt_depth ) / ( 6._r8 * ( 4._r8 + opt_depth ) + opt_depth**2 )
122 0 : radf = max( radinvfrac * qrlw(i,kt) / ( pi(i,kt) - pi(i,kt+1) ) * ( zi(i,kt) - zi(i,kt+1) ), 0._r8 )
123 : ! 2. From 'kt-1' layer and add the contribution from 'kt' layer
124 0 : lwp = ql(i,kt-1) * ( pi(i,kt) - pi(i,kt-1) ) / g
125 0 : opt_depth = 156._r8 * lwp ! Estimated LW optical depth in the CL top layer
126 0 : radinvfrac = opt_depth * ( 4._r8 + opt_depth ) / ( 6._r8 * ( 4._r8 + opt_depth) + opt_depth**2 )
127 0 : radf = radf + max( radinvfrac * qrlw(i,kt-1) / ( pi(i,kt-1) - pi(i,kt) ) * ( zi(i,kt-1) - zi(i,kt) ), 0._r8 )
128 0 : radf = max( radf, 0._r8 ) * chs(i,kt)
129 :
130 : endif
131 :
132 0 : lwp_CL(ncv) = lwp
133 0 : opt_depth_CL(ncv) = opt_depth
134 0 : radinvfrac_CL(ncv) = radinvfrac
135 0 : radf_CL(ncv) = radf
136 : end do ! ncv = 1, ncvfin(i)
137 0 : end subroutine compute_radf
138 :
139 : end module pbl_utils
|