CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
radiation_cloud_overlap.F90
1
3
6 use machine, only : kind_phys
7 implicit none
8 public :: cmp_dcorr_lgth
9
11 module procedure cmp_dcorr_lgth_hogan
12 module procedure cmp_dcorr_lgth_oreopoulos
13 end interface
14
15contains
16
17 ! ######################################################################################
18 ! Hogan et al. (2010)
19 ! "Effect of improving representation of horizontal and vertical cloud structure on the
20 ! Earth's global radiation budget. Part I: Review and parametrization"
21 ! https://rmets.onlinelibrary.wiley.com/doi/full/10.1002/qj.647
22 ! ######################################################################################
24 subroutine cmp_dcorr_lgth_hogan(nCol, lat, con_pi, dcorr_lgth)
25 ! Inputs
26 integer, intent(in) :: &
27 nCol ! Number of horizontal grid-points
28 real(kind_phys), intent(in) :: &
29 con_pi ! Physical constant: Pi
30 real(kind_phys), dimension(:), intent(in) :: &
31 lat ! Latitude
32 ! Outputs
33 real(kind_phys), dimension(:),intent(out) :: &
34 dcorr_lgth ! Decorrelation length
35
36 ! Local variables
37 integer :: iCol
38
39 ! Parameters
40 real(kind_phys),parameter :: min_dcorr = 0.6 ! (see section 2.3)
41
42 do icol =1,ncol
43 dcorr_lgth(icol) = max(min_dcorr, 2.78-4.6*abs(lat(icol)/con_pi)) ! (eq. 13)
44 enddo
45
46 end subroutine cmp_dcorr_lgth_hogan
47 ! ######################################################################################
48 ! Oreopoulos et al. (2012)
49 ! "Radiative impacts of cloud heterogeneity and overlap in an
50 ! atmospheric General Circulation Model"
51 ! 10.5194/acp-12-9097-2012
52 ! ######################################################################################
54 subroutine cmp_dcorr_lgth_oreopoulos(nCol, lat, juldat, yearlength, dcorr_lgth)
55 ! Inputs
56 integer, intent(in) :: &
57 nCol, & ! Number of horizontal grid-points
58 yearlength ! Number of days in year
59
60 real(kind_phys), intent(in) :: &
61 juldat ! Julian date
62 real(kind_phys), dimension(:), intent(in) :: &
63 lat ! Latitude
64
65 ! Outputs
66 real(kind_phys), dimension(:),intent(out) :: &
67 dcorr_lgth ! Decorrelation length (km)
68
69 ! Parameters for the Gaussian fits per Eqs. (10) and (11) (See Table 1)
70 real(kind_phys), parameter :: & !
71 am1 = 1.4315_kind_phys, & !
72 am2 = 2.1219_kind_phys, & !
73 am4 = -25.584_kind_phys, & !
74 amr = 7.0_kind_phys !
75
76 ! Local variables
77 integer :: iCol
78 real(kind_phys) :: am3
79
80 do icol = 1, ncol
81 if (juldat .gt. 181._kind_phys) then
82 am3 = -4._kind_phys * amr * (juldat - 272._kind_phys) / yearlength ! (eq. 11a)
83 else
84 am3 = 4._kind_phys * amr * (juldat - 91._kind_phys) / yearlength ! (eq. 11b)
85 endif
86 dcorr_lgth(icol) = am1 + am2 * exp( -(lat(icol) - am3)**2 / am4**2) ! (eq. 10)
87 enddo
88
89 end subroutine cmp_dcorr_lgth_oreopoulos
90
92 subroutine get_alpha_exper(nCol, nLay, iovr, iovr_exprand, dzlay, &
93 dcorr_lgth, cld_frac, alpha)
94
95 ! Inputs
96 integer, intent(in) :: &
97 nCol, & ! Number of horizontal grid points
98 nLay ! Number of vertical grid points
99 integer, intent(in) :: &
100 iovr, &
101 iovr_exprand
102 real(kind_phys), dimension(:), intent(in) :: &
103 dcorr_lgth ! Decorrelation length (km)
104 real(kind_phys), dimension(:,:), intent(in) :: &
105 dzlay !
106 real(kind_phys), dimension(:,:), intent(in) :: &
107 cld_frac
108
109 ! Outputs
110 real(kind_phys), dimension(:,:) :: &
111 alpha ! Cloud overlap parameter
112
113 ! Local variables
114 integer :: iCol,iLay
115
116 do icol = 1, ncol
117 alpha(icol,1) = 0.0d0
118 do ilay = 2, nlay
119 alpha(icol,ilay) = exp( -(dzlay(icol,ilay)) / dcorr_lgth(icol))
120 enddo
121 enddo
122
123 ! Revise alpha for exponential-random cloud overlap
124 ! Decorrelate layers when a clear layer follows a cloudy layer to enforce
125 ! random correlation between non-adjacent blocks of cloudy layers
126 if (iovr == iovr_exprand) then
127 do ilay = 2, nlay
128 do icol = 1, ncol
129 if (cld_frac(icol,ilay) == 0.0 .and. cld_frac(icol,ilay-1) > 0.0) then
130 alpha(icol,ilay) = 0.0
131 endif
132 enddo
133 enddo
134 endif
135
136 return
137
138 end subroutine get_alpha_exper
This module contains the calculation of cloud overlap parameters for both RRTMG and RRTMGP.