1 | ! radiation_mcica_lw.F90 - Monte-Carlo Independent Column Approximation longtwave solver |
---|
2 | ! |
---|
3 | ! (C) Copyright 2015- ECMWF. |
---|
4 | ! |
---|
5 | ! This software is licensed under the terms of the Apache Licence Version 2.0 |
---|
6 | ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0. |
---|
7 | ! |
---|
8 | ! In applying this licence, ECMWF does not waive the privileges and immunities |
---|
9 | ! granted to it by virtue of its status as an intergovernmental organisation |
---|
10 | ! nor does it submit to any jurisdiction. |
---|
11 | ! |
---|
12 | ! Author: Robin Hogan |
---|
13 | ! Email: r.j.hogan@ecmwf.int |
---|
14 | ! |
---|
15 | ! Modifications |
---|
16 | ! 2017-04-11 R. Hogan Receive emission/albedo rather than planck/emissivity |
---|
17 | ! 2017-04-22 R. Hogan Store surface fluxes at all g-points |
---|
18 | ! 2017-07-12 R. Hogan Call fast adding method if only clouds scatter |
---|
19 | ! 2017-10-23 R. Hogan Renamed single-character variables |
---|
20 | |
---|
21 | #include "ecrad_config.h" |
---|
22 | |
---|
23 | module radiation_mcica_lw |
---|
24 | |
---|
25 | public |
---|
26 | |
---|
27 | contains |
---|
28 | |
---|
29 | !--------------------------------------------------------------------- |
---|
30 | ! Longwave Monte Carlo Independent Column Approximation |
---|
31 | ! (McICA). This implementation performs a clear-sky and a cloudy-sky |
---|
32 | ! calculation, and then weights the two to get the all-sky fluxes |
---|
33 | ! according to the total cloud cover. This method reduces noise for |
---|
34 | ! low cloud cover situations, and exploits the clear-sky |
---|
35 | ! calculations that are usually performed for diagnostic purposes |
---|
36 | ! simultaneously. The cloud generator has been carefully written |
---|
37 | ! such that the stochastic cloud field satisfies the prescribed |
---|
38 | ! overlap parameter accounting for this weighting. |
---|
39 | subroutine solver_mcica_lw(nlev,istartcol,iendcol, & |
---|
40 | & config, single_level, cloud, & |
---|
41 | & od, ssa, g, od_cloud, ssa_cloud, g_cloud, planck_hl, & |
---|
42 | & emission, albedo, & |
---|
43 | & flux) |
---|
44 | |
---|
45 | use parkind1, only : jprb |
---|
46 | use yomhook, only : lhook, dr_hook, jphook |
---|
47 | |
---|
48 | use radiation_io, only : nulerr, radiation_abort |
---|
49 | use radiation_config, only : config_type |
---|
50 | use radiation_single_level, only : single_level_type |
---|
51 | use radiation_cloud, only : cloud_type |
---|
52 | use radiation_flux, only : flux_type |
---|
53 | use radiation_two_stream, only : calc_ref_trans_lw, & |
---|
54 | & calc_no_scattering_transmittance_lw |
---|
55 | use radiation_adding_ica_lw, only : adding_ica_lw, fast_adding_ica_lw, & |
---|
56 | & calc_fluxes_no_scattering_lw |
---|
57 | use radiation_lw_derivatives, only : calc_lw_derivatives_ica, modify_lw_derivatives_ica |
---|
58 | use radiation_cloud_generator, only: cloud_generator |
---|
59 | |
---|
60 | implicit none |
---|
61 | |
---|
62 | ! Inputs |
---|
63 | integer, intent(in) :: nlev ! number of model levels |
---|
64 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
65 | type(config_type), intent(in) :: config |
---|
66 | type(single_level_type), intent(in) :: single_level |
---|
67 | type(cloud_type), intent(in) :: cloud |
---|
68 | |
---|
69 | ! Gas and aerosol optical depth, single-scattering albedo and |
---|
70 | ! asymmetry factor at each longwave g-point |
---|
71 | real(jprb), intent(in), dimension(config%n_g_lw, nlev, istartcol:iendcol) :: & |
---|
72 | & od |
---|
73 | real(jprb), intent(in), dimension(config%n_g_lw_if_scattering, nlev, istartcol:iendcol) :: & |
---|
74 | & ssa, g |
---|
75 | |
---|
76 | ! Cloud and precipitation optical depth, single-scattering albedo and |
---|
77 | ! asymmetry factor in each longwave band |
---|
78 | real(jprb), intent(in), dimension(config%n_bands_lw,nlev,istartcol:iendcol) :: & |
---|
79 | & od_cloud |
---|
80 | real(jprb), intent(in), dimension(config%n_bands_lw_if_scattering, & |
---|
81 | & nlev,istartcol:iendcol) :: ssa_cloud, g_cloud |
---|
82 | |
---|
83 | ! Planck function at each half-level and the surface |
---|
84 | real(jprb), intent(in), dimension(config%n_g_lw,nlev+1,istartcol:iendcol) :: & |
---|
85 | & planck_hl |
---|
86 | |
---|
87 | ! Emission (Planck*emissivity) and albedo (1-emissivity) at the |
---|
88 | ! surface at each longwave g-point |
---|
89 | real(jprb), intent(in), dimension(config%n_g_lw, istartcol:iendcol) :: emission, albedo |
---|
90 | |
---|
91 | ! Output |
---|
92 | type(flux_type), intent(inout):: flux |
---|
93 | |
---|
94 | ! Local variables |
---|
95 | |
---|
96 | ! Diffuse reflectance and transmittance for each layer in clear |
---|
97 | ! and all skies |
---|
98 | real(jprb), dimension(config%n_g_lw, nlev) :: ref_clear, trans_clear, reflectance, transmittance |
---|
99 | |
---|
100 | ! Emission by a layer into the upwelling or downwelling diffuse |
---|
101 | ! streams, in clear and all skies |
---|
102 | real(jprb), dimension(config%n_g_lw, nlev) :: source_up_clear, source_dn_clear, source_up, source_dn |
---|
103 | |
---|
104 | ! Fluxes per g point |
---|
105 | real(jprb), dimension(config%n_g_lw, nlev+1) :: flux_up, flux_dn |
---|
106 | real(jprb), dimension(config%n_g_lw, nlev+1) :: flux_up_clear, flux_dn_clear |
---|
107 | |
---|
108 | ! Combined gas+aerosol+cloud optical depth, single scattering |
---|
109 | ! albedo and asymmetry factor |
---|
110 | real(jprb), dimension(config%n_g_lw) :: od_total, ssa_total, g_total |
---|
111 | |
---|
112 | ! Combined scattering optical depth |
---|
113 | real(jprb) :: scat_od, scat_od_total(config%n_g_lw) |
---|
114 | |
---|
115 | ! Optical depth scaling from the cloud generator, zero indicating |
---|
116 | ! clear skies |
---|
117 | real(jprb), dimension(config%n_g_lw,nlev) :: od_scaling |
---|
118 | |
---|
119 | ! Modified optical depth after McICA scaling to represent cloud |
---|
120 | ! inhomogeneity |
---|
121 | real(jprb), dimension(config%n_g_lw) :: od_cloud_new |
---|
122 | |
---|
123 | ! Total cloud cover output from the cloud generator |
---|
124 | real(jprb) :: total_cloud_cover |
---|
125 | |
---|
126 | ! Identify clear-sky layers |
---|
127 | logical :: is_clear_sky_layer(nlev) |
---|
128 | |
---|
129 | ! Temporary storage for more efficient summation |
---|
130 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
131 | real(jprb), dimension(nlev+1,2) :: sum_aux |
---|
132 | #else |
---|
133 | real(jprb) :: sum_up, sum_dn |
---|
134 | #endif |
---|
135 | |
---|
136 | ! Index of the highest cloudy layer |
---|
137 | integer :: i_cloud_top |
---|
138 | |
---|
139 | ! Number of g points |
---|
140 | integer :: ng |
---|
141 | |
---|
142 | ! Loop indices for level, column and g point |
---|
143 | integer :: jlev, jcol, jg |
---|
144 | |
---|
145 | real(jphook) :: hook_handle |
---|
146 | |
---|
147 | if (lhook) call dr_hook('radiation_mcica_lw:solver_mcica_lw',0,hook_handle) |
---|
148 | |
---|
149 | if (.not. config%do_clear) then |
---|
150 | write(nulerr,'(a)') '*** Error: longwave McICA requires clear-sky calculation to be performed' |
---|
151 | call radiation_abort() |
---|
152 | end if |
---|
153 | |
---|
154 | ng = config%n_g_lw |
---|
155 | |
---|
156 | ! Loop through columns |
---|
157 | do jcol = istartcol,iendcol |
---|
158 | |
---|
159 | ! Clear-sky calculation |
---|
160 | if (config%do_lw_aerosol_scattering) then |
---|
161 | ! Scattering case: first compute clear-sky reflectance, |
---|
162 | ! transmittance etc at each model level |
---|
163 | call calc_ref_trans_lw(ng*nlev, & |
---|
164 | & od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), & |
---|
165 | & planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1,jcol), & |
---|
166 | & ref_clear, trans_clear, & |
---|
167 | & source_up_clear, source_dn_clear) |
---|
168 | ! Then use adding method to compute fluxes |
---|
169 | call adding_ica_lw(ng, nlev, & |
---|
170 | & ref_clear, trans_clear, source_up_clear, source_dn_clear, & |
---|
171 | & emission(:,jcol), albedo(:,jcol), & |
---|
172 | & flux_up_clear, flux_dn_clear) |
---|
173 | else |
---|
174 | ! Non-scattering case: use simpler functions for |
---|
175 | ! transmission and emission |
---|
176 | call calc_no_scattering_transmittance_lw(ng*nlev, od(:,:,jcol), & |
---|
177 | & planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1, jcol), & |
---|
178 | & trans_clear, source_up_clear, source_dn_clear) |
---|
179 | ! Ensure that clear-sky reflectance is zero since it may be |
---|
180 | ! used in cloudy-sky case |
---|
181 | ref_clear = 0.0_jprb |
---|
182 | ! Simpler down-then-up method to compute fluxes |
---|
183 | call calc_fluxes_no_scattering_lw(ng, nlev, & |
---|
184 | & trans_clear, source_up_clear, source_dn_clear, & |
---|
185 | & emission(:,jcol), albedo(:,jcol), & |
---|
186 | & flux_up_clear, flux_dn_clear) |
---|
187 | end if |
---|
188 | |
---|
189 | ! Sum over g-points to compute broadband fluxes |
---|
190 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
191 | sum_aux(:,:) = 0.0_jprb |
---|
192 | do jg = 1,ng |
---|
193 | do jlev = 1,nlev+1 |
---|
194 | sum_aux(jlev,1) = sum_aux(jlev,1) + flux_up_clear(jg,jlev) |
---|
195 | sum_aux(jlev,2) = sum_aux(jlev,2) + flux_dn_clear(jg,jlev) |
---|
196 | end do |
---|
197 | end do |
---|
198 | flux%lw_up_clear(jcol,:) = sum_aux(:,1) |
---|
199 | flux%lw_dn_clear(jcol,:) = sum_aux(:,2) |
---|
200 | #else |
---|
201 | do jlev = 1,nlev+1 |
---|
202 | sum_up = 0.0_jprb |
---|
203 | sum_dn = 0.0_jprb |
---|
204 | !$omp simd reduction(+:sum_up, sum_dn) |
---|
205 | do jg = 1,ng |
---|
206 | sum_up = sum_up + flux_up_clear(jg,jlev) |
---|
207 | sum_dn = sum_dn + flux_dn_clear(jg,jlev) |
---|
208 | end do |
---|
209 | flux%lw_up_clear(jcol,jlev) = sum_up |
---|
210 | flux%lw_dn_clear(jcol,jlev) = sum_dn |
---|
211 | end do |
---|
212 | #endif |
---|
213 | |
---|
214 | ! Store surface spectral downwelling fluxes |
---|
215 | flux%lw_dn_surf_clear_g(:,jcol) = flux_dn_clear(:,nlev+1) |
---|
216 | |
---|
217 | ! Do cloudy-sky calculation; add a prime number to the seed in |
---|
218 | ! the longwave |
---|
219 | call cloud_generator(ng, nlev, config%i_overlap_scheme, & |
---|
220 | & single_level%iseed(jcol) + 997, & |
---|
221 | & config%cloud_fraction_threshold, & |
---|
222 | & cloud%fraction(jcol,:), cloud%overlap_param(jcol,:), & |
---|
223 | & config%cloud_inhom_decorr_scaling, cloud%fractional_std(jcol,:), & |
---|
224 | & config%pdf_sampler, od_scaling, total_cloud_cover, & |
---|
225 | & use_beta_overlap=config%use_beta_overlap, & |
---|
226 | & use_vectorizable_generator=config%use_vectorizable_generator) |
---|
227 | |
---|
228 | ! Store total cloud cover |
---|
229 | flux%cloud_cover_lw(jcol) = total_cloud_cover |
---|
230 | |
---|
231 | if (total_cloud_cover >= config%cloud_fraction_threshold) then |
---|
232 | ! Total-sky calculation |
---|
233 | |
---|
234 | is_clear_sky_layer = .true. |
---|
235 | i_cloud_top = nlev+1 |
---|
236 | do jlev = 1,nlev |
---|
237 | ! Compute combined gas+aerosol+cloud optical properties |
---|
238 | if (cloud%fraction(jcol,jlev) >= config%cloud_fraction_threshold) then |
---|
239 | is_clear_sky_layer(jlev) = .false. |
---|
240 | ! Get index to the first cloudy layer from the top |
---|
241 | if (i_cloud_top > jlev) then |
---|
242 | i_cloud_top = jlev |
---|
243 | end if |
---|
244 | |
---|
245 | do jg = 1,ng |
---|
246 | od_cloud_new(jg) = od_scaling(jg,jlev) & |
---|
247 | & * od_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) |
---|
248 | od_total(jg) = od(jg,jlev,jcol) + od_cloud_new(jg) |
---|
249 | ssa_total(jg) = 0.0_jprb |
---|
250 | g_total(jg) = 0.0_jprb |
---|
251 | end do |
---|
252 | |
---|
253 | if (config%do_lw_cloud_scattering) then |
---|
254 | ! Scattering case: calculate reflectance and |
---|
255 | ! transmittance at each model level |
---|
256 | |
---|
257 | if (config%do_lw_aerosol_scattering) then |
---|
258 | ! In single precision we need to protect against the |
---|
259 | ! case that od_total > 0.0 and ssa_total > 0.0 but |
---|
260 | ! od_total*ssa_total == 0 due to underflow |
---|
261 | do jg = 1,ng |
---|
262 | if (od_total(jg) > 0.0_jprb) then |
---|
263 | scat_od_total(jg) = ssa(jg,jlev,jcol)*od(jg,jlev,jcol) & |
---|
264 | & + ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
265 | & * od_cloud_new(jg) |
---|
266 | ssa_total(jg) = scat_od_total(jg) / od_total(jg) |
---|
267 | |
---|
268 | if (scat_od_total(jg) > 0.0_jprb) then |
---|
269 | g_total(jg) = (g(jg,jlev,jcol)*ssa(jg,jlev,jcol)*od(jg,jlev,jcol) & |
---|
270 | & + g_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
271 | & * ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
272 | & * od_cloud_new(jg)) & |
---|
273 | & / scat_od_total(jg) |
---|
274 | end if |
---|
275 | end if |
---|
276 | end do |
---|
277 | |
---|
278 | else |
---|
279 | |
---|
280 | do jg = 1,ng |
---|
281 | if (od_total(jg) > 0.0_jprb) then |
---|
282 | scat_od = ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
283 | & * od_cloud_new(jg) |
---|
284 | ssa_total(jg) = scat_od / od_total(jg) |
---|
285 | if (scat_od > 0.0_jprb) then |
---|
286 | g_total(jg) = g_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
287 | & * ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) & |
---|
288 | & * od_cloud_new(jg) / scat_od |
---|
289 | end if |
---|
290 | end if |
---|
291 | end do |
---|
292 | |
---|
293 | end if |
---|
294 | |
---|
295 | ! Compute cloudy-sky reflectance, transmittance etc at |
---|
296 | ! each model level |
---|
297 | call calc_ref_trans_lw(ng, & |
---|
298 | & od_total, ssa_total, g_total, & |
---|
299 | & planck_hl(:,jlev,jcol), planck_hl(:,jlev+1,jcol), & |
---|
300 | & reflectance(:,jlev), transmittance(:,jlev), & |
---|
301 | & source_up(:,jlev), source_dn(:,jlev)) |
---|
302 | else |
---|
303 | ! No-scattering case: use simpler functions for |
---|
304 | ! transmission and emission |
---|
305 | call calc_no_scattering_transmittance_lw(ng, od_total, & |
---|
306 | & planck_hl(:,jlev,jcol), planck_hl(:,jlev+1, jcol), & |
---|
307 | & transmittance(:,jlev), source_up(:,jlev), source_dn(:,jlev)) |
---|
308 | end if |
---|
309 | |
---|
310 | else |
---|
311 | ! Clear-sky layer: copy over clear-sky values |
---|
312 | do jg = 1,ng |
---|
313 | reflectance(jg,jlev) = ref_clear(jg,jlev) |
---|
314 | transmittance(jg,jlev) = trans_clear(jg,jlev) |
---|
315 | source_up(jg,jlev) = source_up_clear(jg,jlev) |
---|
316 | source_dn(jg,jlev) = source_dn_clear(jg,jlev) |
---|
317 | end do |
---|
318 | end if |
---|
319 | end do |
---|
320 | |
---|
321 | if (config%do_lw_aerosol_scattering) then |
---|
322 | ! Use adding method to compute fluxes for an overcast sky, |
---|
323 | ! allowing for scattering in all layers |
---|
324 | call adding_ica_lw(ng, nlev, reflectance, transmittance, source_up, source_dn, & |
---|
325 | & emission(:,jcol), albedo(:,jcol), & |
---|
326 | & flux_up, flux_dn) |
---|
327 | else if (config%do_lw_cloud_scattering) then |
---|
328 | ! Use adding method to compute fluxes but optimize for the |
---|
329 | ! presence of clear-sky layers |
---|
330 | call fast_adding_ica_lw(ng, nlev, reflectance, transmittance, source_up, source_dn, & |
---|
331 | & emission(:,jcol), albedo(:,jcol), & |
---|
332 | & is_clear_sky_layer, i_cloud_top, flux_dn_clear, & |
---|
333 | & flux_up, flux_dn) |
---|
334 | else |
---|
335 | ! Simpler down-then-up method to compute fluxes |
---|
336 | call calc_fluxes_no_scattering_lw(ng, nlev, & |
---|
337 | & transmittance, source_up, source_dn, emission(:,jcol), albedo(:,jcol), & |
---|
338 | & flux_up, flux_dn) |
---|
339 | end if |
---|
340 | |
---|
341 | ! Store overcast broadband fluxes |
---|
342 | #ifdef DWD_REDUCTION_OPTIMIZATIONS |
---|
343 | sum_aux(:,:) = 0._jprb |
---|
344 | do jg = 1, ng |
---|
345 | do jlev = 1, nlev+1 |
---|
346 | sum_aux(jlev,1) = sum_aux(jlev,1) + flux_up(jg,jlev) |
---|
347 | sum_aux(jlev,2) = sum_aux(jlev,2) + flux_dn(jg,jlev) |
---|
348 | end do |
---|
349 | end do |
---|
350 | flux%lw_up(jcol,:) = sum_aux(:,1) |
---|
351 | flux%lw_dn(jcol,:) = sum_aux(:,2) |
---|
352 | #else |
---|
353 | do jlev = 1,nlev+1 |
---|
354 | sum_up = 0.0_jprb |
---|
355 | sum_dn = 0.0_jprb |
---|
356 | !$omp simd reduction(+:sum_up, sum_dn) |
---|
357 | do jg = 1,ng |
---|
358 | sum_up = sum_up + flux_up(jg,jlev) |
---|
359 | sum_dn = sum_dn + flux_dn(jg,jlev) |
---|
360 | end do |
---|
361 | flux%lw_up(jcol,jlev) = sum_up |
---|
362 | flux%lw_dn(jcol,jlev) = sum_dn |
---|
363 | end do |
---|
364 | #endif |
---|
365 | |
---|
366 | ! Cloudy flux profiles currently assume completely overcast |
---|
367 | ! skies; perform weighted average with clear-sky profile |
---|
368 | do jlev = 1,nlev+1 |
---|
369 | flux%lw_up(jcol,jlev) = total_cloud_cover *flux%lw_up(jcol,jlev) & |
---|
370 | & + (1.0_jprb - total_cloud_cover)*flux%lw_up_clear(jcol,jlev) |
---|
371 | flux%lw_dn(jcol,jlev) = total_cloud_cover *flux%lw_dn(jcol,jlev) & |
---|
372 | & + (1.0_jprb - total_cloud_cover)*flux%lw_dn_clear(jcol,jlev) |
---|
373 | end do |
---|
374 | ! Store surface spectral downwelling fluxes |
---|
375 | flux%lw_dn_surf_g(:,jcol) = total_cloud_cover*flux_dn(:,nlev+1) & |
---|
376 | & + (1.0_jprb - total_cloud_cover)*flux%lw_dn_surf_clear_g(:,jcol) |
---|
377 | |
---|
378 | ! Compute the longwave derivatives needed by Hogan and Bozzo |
---|
379 | ! (2015) approximate radiation update scheme |
---|
380 | if (config%do_lw_derivatives) then |
---|
381 | call calc_lw_derivatives_ica(ng, nlev, jcol, transmittance, flux_up(:,nlev+1), & |
---|
382 | & flux%lw_derivatives) |
---|
383 | if (total_cloud_cover < 1.0_jprb - config%cloud_fraction_threshold) then |
---|
384 | ! Modify the existing derivative with the contribution from the clear sky |
---|
385 | call modify_lw_derivatives_ica(ng, nlev, jcol, trans_clear, flux_up_clear(:,nlev+1), & |
---|
386 | & 1.0_jprb-total_cloud_cover, flux%lw_derivatives) |
---|
387 | end if |
---|
388 | end if |
---|
389 | |
---|
390 | else |
---|
391 | ! No cloud in profile and clear-sky fluxes already |
---|
392 | ! calculated: copy them over |
---|
393 | do jlev = 1,nlev+1 |
---|
394 | flux%lw_up(jcol,jlev) = flux%lw_up_clear(jcol,jlev) |
---|
395 | flux%lw_dn(jcol,jlev) = flux%lw_dn_clear(jcol,jlev) |
---|
396 | end do |
---|
397 | flux%lw_dn_surf_g(:,jcol) = flux%lw_dn_surf_clear_g(:,jcol) |
---|
398 | if (config%do_lw_derivatives) then |
---|
399 | call calc_lw_derivatives_ica(ng, nlev, jcol, trans_clear, flux_up_clear(:,nlev+1), & |
---|
400 | & flux%lw_derivatives) |
---|
401 | |
---|
402 | end if |
---|
403 | end if ! Cloud is present in profile |
---|
404 | end do |
---|
405 | |
---|
406 | if (lhook) call dr_hook('radiation_mcica_lw:solver_mcica_lw',1,hook_handle) |
---|
407 | |
---|
408 | end subroutine solver_mcica_lw |
---|
409 | |
---|
410 | end module radiation_mcica_lw |
---|