36 function draw_samples(cloud_mask,do_twostream,clouds,clouds_sampled)
result(error_msg)
38 logical,
dimension(:,:,:),
intent(in ) :: cloud_mask
39 logical,
intent(in ) :: do_twostream
40 class(ty_optical_props_arry),
intent(in ) :: clouds
43 class(ty_optical_props_arry),
intent(inout) :: clouds_sampled
44 character(len=128) :: error_msg
47 integer :: ncol,nlay,nbnd,ngpt
53 ncol = clouds%get_ncol()
54 nlay = clouds%get_nlay()
55 nbnd = clouds%get_nband()
56 ngpt = clouds_sampled%get_ngpt()
59 call apply_cloud_mask(ncol,nlay,nbnd,ngpt,clouds_sampled%get_band_lims_gpoint(),cloud_mask,clouds%tau,clouds_sampled%tau)
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 )
72 error_msg =
"draw_samples: by-band and sampled cloud properties need to be the same variable type"
78 subroutine sampled_mask(randoms, cloud_frac, cloud_mask, overlap_param, randoms2)
80 real(wp),
dimension(:,:,:),
intent(in ) :: randoms
81 real(wp),
dimension(:,:),
intent(in ) :: cloud_frac
84 logical,
dimension(:,:,:),
intent(out) :: cloud_mask
87 real(wp),
dimension(:,:),
intent(in ),
optional :: overlap_param
88 real(wp),
dimension(:,:,:),
intent(in ),
optional :: randoms2
91 integer :: ncol, nlay, ngpt, icol, ilay, igpt
92 integer :: cloud_lay_fst, cloud_lay_lst
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
101 ncol =
size(randoms, 3)
102 nlay =
size(randoms, 2)
103 ngpt =
size(randoms, 1)
106 if (
present(overlap_param)) l_use_overlap_param = .true.
109 if (
present(randoms2)) l_use_second_rng = .true.
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.
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.
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
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))
136 cloud_mask(icol,ilay,1:ngpt) = .false.
144 if (l_use_overlap_param)
then
145 if(cloud_mask_layer(ilay))
then
146 if(cloud_mask_layer(ilay-1))
then
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
153 local_rands(1:ngpt) = randoms(1:ngpt,ilay,icol)
155 cloud_mask(icol,ilay,1:ngpt) = local_rands(1:ngpt) > (1._wp - cloud_frac(icol,ilay))
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))
167 cloud_mask(icol,ilay,1:ngpt) = randoms(1:ngpt,ilay,icol) > (1._wp - cloud_frac(icol,ilay))
173 cloud_mask(icol,cloud_lay_lst+1:nlay, 1:ngpt) = .false.
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
187 integer :: icol,ilay,ibnd,igpt
190 do igpt = band_lims_gpt(1,ibnd), band_lims_gpt(2,ibnd)
192 sampled_field(1:ncol,ilay,igpt) = merge(input_field(1:ncol,ilay,ibnd), 0._wp, cloud_mask(1:ncol,ilay,igpt))