CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
rrtmgp_sampling.F90
1
3!
4! Contacts: Robert Pincus and Eli Mlawer
5! email: rrtmgp@aer.com
6!
7! Copyright 2015-2019, Atmospheric and Environmental Research and
8! Regents of the University of Colorado. All right reserved.
9!
10! Use and duplication is permitted under the terms of the
11! BSD 3-clause license, see http://opensource.org/licenses/BSD-3-Clause
12! -------------------------------------------------------------------------------------------------
13!
16! Cloud optical properties, defined by band and assumed homogenous within each cell (column/layer),
17! are randomly sampled to preserve the mean cloud fraction and one of several possible overlap assumptions
18! Users supply random numbers with order ngpt,nlay,ncol
19! These are only accessed if cloud_fraction(icol,ilay) > 0 so many values don't need to be filled in
20!
21! Adapted by Dustin Swales on 8/11/2020 for use in UFS (NOAA-PSL/CU-CIRES)
22!
23! -------------------------------------------------------------------------------------------------
25 use mo_rte_kind, only: wp, wl
26 use mo_optical_props, only: ty_optical_props_arry, &
27 ty_optical_props_1scl, &
28 ty_optical_props_2str, &
29 ty_optical_props_nstr
30 implicit none
31 private
32 public :: draw_samples, sampled_mask
33contains
36 function draw_samples(cloud_mask,do_twostream,clouds,clouds_sampled) result(error_msg)
37 ! Inputs
38 logical, dimension(:,:,:), intent(in ) :: cloud_mask ! Dimensions ncol,nlay,ngpt
39 logical, intent(in ) :: do_twostream ! Do two-stream?
40 class(ty_optical_props_arry), intent(in ) :: clouds ! Defined by band
41
42 ! Outputs
43 class(ty_optical_props_arry), intent(inout) :: clouds_sampled ! Defined by g-point
44 character(len=128) :: error_msg
45
46 ! Local variables
47 integer :: ncol,nlay,nbnd,ngpt
48 integer :: imom
49
50 error_msg = ""
51
52 ! Array extents
53 ncol = clouds%get_ncol()
54 nlay = clouds%get_nlay()
55 nbnd = clouds%get_nband()
56 ngpt = clouds_sampled%get_ngpt()
57
58 ! Optical depth assignment works for 1scl, 2str (also nstr)
59 call apply_cloud_mask(ncol,nlay,nbnd,ngpt,clouds_sampled%get_band_lims_gpoint(),cloud_mask,clouds%tau,clouds_sampled%tau)
60 !
61 ! For 2-stream
62 !
63 select type(clouds)
64 type is (ty_optical_props_2str)
65 select type(clouds_sampled)
66 type is (ty_optical_props_2str)
67 if (do_twostream) then
68 call apply_cloud_mask(ncol,nlay,nbnd,ngpt,clouds_sampled%get_band_lims_gpoint(),cloud_mask,clouds%ssa,clouds_sampled%ssa)
69 call apply_cloud_mask(ncol,nlay,nbnd,ngpt,clouds_sampled%get_band_lims_gpoint(),cloud_mask,clouds%g, clouds_sampled%g )
70 endif
71 class default
72 error_msg = "draw_samples: by-band and sampled cloud properties need to be the same variable type"
73 end select
74 end select
75 end function draw_samples
76
78 subroutine sampled_mask(randoms, cloud_frac, cloud_mask, overlap_param, randoms2)
79 ! Inputs
80 real(wp), dimension(:,:,:), intent(in ) :: randoms ! ngpt,nlay,ncol
81 real(wp), dimension(:,:), intent(in ) :: cloud_frac ! ncol,nlay
82
83 ! Outputs
84 logical, dimension(:,:,:), intent(out) :: cloud_mask ! ncol,nlay,ngpt
85
86 ! Inputs (optional)
87 real(wp), dimension(:,:), intent(in ), optional :: overlap_param ! ncol,nlay-1
88 real(wp), dimension(:,:,:), intent(in ), optional :: randoms2 ! ngpt,nlay,ncol
89
90 ! Local variables
91 integer :: ncol, nlay, ngpt, icol, ilay, igpt
92 integer :: cloud_lay_fst, cloud_lay_lst
93 real(wp) :: rho
94 real(wp), dimension(size(randoms,1)) :: local_rands
95 logical, dimension(size(randoms,2)) :: cloud_mask_layer
96 logical :: l_use_overlap_param = .false.
97 logical :: l_use_second_rng = .false.
98 character(len=128) :: error_msg
99
100 ! Array dimensions
101 ncol = size(randoms, 3)
102 nlay = size(randoms, 2)
103 ngpt = size(randoms, 1)
104
105 ! Using cloud-overlap parameter (alpha)?
106 if (present(overlap_param)) l_use_overlap_param = .true.
107
108 ! Using a second RNG?
109 if (present(randoms2)) l_use_second_rng = .true.
110
111 ! Construct the cloud mask for each column
112 do icol = 1, ncol
113 cloud_mask_layer(1:nlay) = cloud_frac(icol,1:nlay) > 0._wp
114 if(.not. any(cloud_mask_layer)) then
115 cloud_mask(icol,1:nlay,1:ngpt) = .false.
116 cycle
117 end if
118 cloud_lay_fst = findloc(cloud_mask_layer, .true., dim=1)
119 cloud_lay_lst = findloc(cloud_mask_layer, .true., dim=1, back = .true.)
120 cloud_mask(icol,1:cloud_lay_fst,1:ngpt) = .false.
121
122 ilay = cloud_lay_fst
123 local_rands(1:ngpt) = randoms(1:ngpt,cloud_lay_fst,icol)
124 cloud_mask(icol,ilay,1:ngpt) = local_rands(1:ngpt) > (1._wp - cloud_frac(icol,ilay))
125 do ilay = cloud_lay_fst+1, cloud_lay_lst
126 ! ################################################################################
127 ! Max-random overlap
128 ! new random deviates if the adjacent layer isn't cloudy
129 ! same random deviates if the adjacent layer is cloudy
130 ! ################################################################################
131 if (.not. l_use_overlap_param) then
132 if(cloud_mask_layer(ilay)) then
133 if(.not. cloud_mask_layer(ilay-1)) local_rands(1:ngpt) = randoms(1:ngpt,ilay,icol)
134 cloud_mask(icol,ilay,1:ngpt) = local_rands(1:ngpt) > (1._wp - cloud_frac(icol,ilay))
135 else
136 cloud_mask(icol,ilay,1:ngpt) = .false.
137 end if
138 end if ! END COND: Maximum-random overlap
139 ! ################################################################################
140 ! Exponential-random overlap
141 ! new random deviates if the adjacent layer isn't cloudy
142 ! correlated deviates if the adjacent layer is cloudy
143 ! ################################################################################
144 if (l_use_overlap_param) then
145 if(cloud_mask_layer(ilay)) then
146 if(cloud_mask_layer(ilay-1)) then
147 ! Create random deviates correlated between this layer and the previous layer
148 ! (have to remove mean value before enforcing correlation).
149 rho = overlap_param(icol,ilay-1)
150 local_rands(1:ngpt) = rho*(local_rands(1:ngpt) -0.5_wp) + &
151 sqrt(1._wp-rho*rho)*(randoms(1:ngpt,ilay,icol)-0.5_wp) + 0.5_wp
152 else
153 local_rands(1:ngpt) = randoms(1:ngpt,ilay,icol)
154 end if
155 cloud_mask(icol,ilay,1:ngpt) = local_rands(1:ngpt) > (1._wp - cloud_frac(icol,ilay))
156 endif
157 endif ! END COND: Exponential/Exponential-random overlap
158 ! ################################################################################
159 ! Exponential-decorrelation overlap
160 ! new random deviates if the adjacent layer isn't cloudy
161 ! correlated deviates if the adjacent layer is cloudy and decorrelation-length
162 ! ################################################################################
163 if (l_use_overlap_param .and. l_use_second_rng) then
164 where(randoms2(1:ngpt,ilay,icol) .le. overlap_param(icol,ilay-1))
165 cloud_mask(icol,ilay,1:ngpt) = randoms(1:ngpt,ilay-1,icol) > (1._wp - cloud_frac(icol,ilay))
166 elsewhere
167 cloud_mask(icol,ilay,1:ngpt) = randoms(1:ngpt,ilay,icol) > (1._wp - cloud_frac(icol,ilay))
168 end where
169 endif ! END COND: Exponential decorrelation-length
170 end do ! END LOOP: Layers
171
172 ! Set cloud-mask in layer below clouds to false
173 cloud_mask(icol,cloud_lay_lst+1:nlay, 1:ngpt) = .false.
174 end do ! END LOOP: Columns
175
176 end subroutine sampled_mask
177
180 subroutine apply_cloud_mask(ncol,nlay,nbnd,ngpt,band_lims_gpt,cloud_mask,input_field,sampled_field)
181 integer, intent(in ) :: ncol,nlay,nbnd,ngpt
182 integer, dimension(2,nbnd), intent(in ) :: band_lims_gpt
183 logical, dimension(ncol,nlay,ngpt), intent(in ) :: cloud_mask
184 real(wp), dimension(ncol,nlay,nbnd), intent(in ) :: input_field
185 real(wp), dimension(ncol,nlay,ngpt), intent(out) :: sampled_field
186
187 integer :: icol,ilay,ibnd,igpt
188
189 do ibnd = 1, nbnd
190 do igpt = band_lims_gpt(1,ibnd), band_lims_gpt(2,ibnd)
191 do ilay = 1, nlay
192 sampled_field(1:ncol,ilay,igpt) = merge(input_field(1:ncol,ilay,ibnd), 0._wp, cloud_mask(1:ncol,ilay,igpt))
193 end do
194 end do
195 end do
196 end subroutine apply_cloud_mask
197
198end module rrtmgp_sampling
This module provides a simple implementation of sampling for the Monte Carlo Independent Pixel Approx...