1 | ! radiation_spartacus_sw.F90 - SPARTACUS shortwave solver |
---|
2 | ! |
---|
3 | ! (C) Copyright 2014- 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 albedos at g-points |
---|
17 | ! 2017-04-22 R. Hogan Store surface fluxes at all g-points |
---|
18 | ! 2017-06-30 R. Hogan Reformulate to use total_albedo_direct not total_source |
---|
19 | ! 2017-07-03 R. Hogan Explicit calculation of encroachment |
---|
20 | ! 2017-10-23 R. Hogan Renamed single-character variables |
---|
21 | ! 2018-02-20 R. Hogan Corrected "computed" encroachment |
---|
22 | ! 2018-03-09 R. Hogan Security on computed encroachment, transmittance and reflectance |
---|
23 | ! 2018-08-29 R. Hogan Reformulate horizontal migration distances in step_migrations |
---|
24 | ! 2018-09-02 R. Hogan Bug fix in x_direct in step_migrations |
---|
25 | ! 2018-09-03 R. Hogan Security via min_cloud_effective_size |
---|
26 | ! 2018-09-04 R. Hogan Use encroachment_scaling and read encroachment edge_length from upper layer |
---|
27 | ! 2018-09-13 R. Hogan Added "Fractal" encroachment option |
---|
28 | ! 2018-09-14 R. Hogan Added "Zero" encroachment option |
---|
29 | ! 2018-10-08 R. Hogan Call calc_region_properties |
---|
30 | ! 2018-10-15 R. Hogan Added call to fast_expm_exchange instead of expm |
---|
31 | ! 2019-01-12 R. Hogan Use inv_inhom_effective_size if allocated |
---|
32 | ! 2019-02-10 R. Hogan Renamed "encroachment" to "entrapment" |
---|
33 | |
---|
34 | module radiation_spartacus_sw |
---|
35 | |
---|
36 | public |
---|
37 | |
---|
38 | contains |
---|
39 | |
---|
40 | ! Small routine for scaling cloud optical depth in the cloudy |
---|
41 | ! regions |
---|
42 | #include "radiation_optical_depth_scaling.h" |
---|
43 | |
---|
44 | ! This module contains just one exported subroutine, the shortwave |
---|
45 | ! solver using the Speedy Algorithm for Radiative Transfer through |
---|
46 | ! Cloud Sides (SPARTACUS), which can represent 3D effects using a |
---|
47 | ! matrix form of the two-stream equations. |
---|
48 | ! |
---|
49 | ! Sections: |
---|
50 | ! 1: Prepare general variables and arrays |
---|
51 | ! 2: Prepare column-specific variables and arrays |
---|
52 | ! 3: First loop over layers |
---|
53 | ! 3.1: Layer-specific initialization |
---|
54 | ! 3.2: Compute gamma variables |
---|
55 | ! 3.2a: Clear-sky case |
---|
56 | ! 3.2b: Cloudy case |
---|
57 | ! 3.3: Compute reflection, transmission and sources |
---|
58 | ! 3.3a: g-points with 3D effects |
---|
59 | ! 3.3b: g-points without 3D effects |
---|
60 | ! 4: Compute total albedos |
---|
61 | ! 4.1: Adding method |
---|
62 | ! 4.2: Overlap and entrapment |
---|
63 | ! 5: Compute fluxes |
---|
64 | subroutine solver_spartacus_sw(nlev,istartcol,iendcol, & |
---|
65 | & config, single_level, thermodynamics, cloud, & |
---|
66 | & od, ssa, g, od_cloud, ssa_cloud, g_cloud, & |
---|
67 | & albedo_direct, albedo_diffuse, incoming_sw, & |
---|
68 | & flux) |
---|
69 | |
---|
70 | use parkind1, only : jprb |
---|
71 | use yomhook, only : lhook, dr_hook, jphook |
---|
72 | |
---|
73 | use radiation_io, only : nulout |
---|
74 | use radiation_config, only : config_type, IPdfShapeGamma, & |
---|
75 | & IEntrapmentZero, IEntrapmentEdgeOnly, IEntrapmentExplicit, & |
---|
76 | & IEntrapmentExplicitNonFractal, IEntrapmentMaximum |
---|
77 | use radiation_single_level, only : single_level_type |
---|
78 | use radiation_thermodynamics, only : thermodynamics_type |
---|
79 | use radiation_cloud, only : cloud_type |
---|
80 | use radiation_regions, only : calc_region_properties |
---|
81 | use radiation_overlap, only : calc_overlap_matrices |
---|
82 | use radiation_flux, only : flux_type, & |
---|
83 | & indexed_sum, add_indexed_sum |
---|
84 | use radiation_matrix |
---|
85 | use radiation_two_stream, only : calc_two_stream_gammas_sw, & |
---|
86 | & calc_reflectance_transmittance_sw, calc_frac_scattered_diffuse_sw |
---|
87 | use radiation_constants, only : Pi, GasConstantDryAir, & |
---|
88 | & AccelDueToGravity |
---|
89 | |
---|
90 | implicit none |
---|
91 | |
---|
92 | ! Inputs |
---|
93 | integer, intent(in) :: nlev ! number of model levels |
---|
94 | integer, intent(in) :: istartcol, iendcol ! range of columns to process |
---|
95 | type(config_type), intent(in) :: config |
---|
96 | type(single_level_type), intent(in) :: single_level |
---|
97 | type(thermodynamics_type),intent(in) :: thermodynamics |
---|
98 | type(cloud_type), intent(in) :: cloud |
---|
99 | |
---|
100 | ! Gas and aerosol optical depth, single-scattering albedo and |
---|
101 | ! asymmetry factor at each shortwave g-point |
---|
102 | real(jprb), intent(in), dimension(config%n_g_sw,nlev,istartcol:iendcol) :: & |
---|
103 | & od, ssa, g |
---|
104 | |
---|
105 | ! Cloud and precipitation optical depth, single-scattering albedo and |
---|
106 | ! asymmetry factor in each shortwave band |
---|
107 | real(jprb), intent(in), dimension(config%n_bands_sw,nlev,istartcol:iendcol) :: & |
---|
108 | & od_cloud, ssa_cloud, g_cloud |
---|
109 | |
---|
110 | ! Direct and diffuse surface albedos, and the incoming shortwave |
---|
111 | ! flux into a plane perpendicular to the incoming radiation at |
---|
112 | ! top-of-atmosphere in each of the shortwave g points |
---|
113 | real(jprb), intent(in), dimension(config%n_g_sw,istartcol:iendcol) :: & |
---|
114 | & albedo_direct, albedo_diffuse, incoming_sw |
---|
115 | |
---|
116 | ! Output |
---|
117 | type(flux_type), intent(inout):: flux |
---|
118 | |
---|
119 | integer :: nreg, ng |
---|
120 | integer :: nregactive ! =1 in clear layer, =nreg in a cloudy layer |
---|
121 | integer :: jcol, jlev, jg, jreg, iband, jreg2,jreg3 |
---|
122 | #ifdef EXPLICIT_EDGE_ENTRAPMENT |
---|
123 | integer :: jreg4 |
---|
124 | #endif |
---|
125 | integer :: ng3D ! Number of g-points with small enough gas optical |
---|
126 | ! depth that 3D effects need to be represented |
---|
127 | |
---|
128 | ! Ratio of gas constant for dry air to acceleration due to gravity |
---|
129 | real(jprb), parameter :: R_over_g = GasConstantDryAir / AccelDueToGravity |
---|
130 | |
---|
131 | real(jprb) :: mu0, one_over_mu0, tan_sza, transfer_scaling |
---|
132 | |
---|
133 | ! The tangent of the effective zenith angle for diffuse radiation |
---|
134 | ! that is appropriate for 3D transport |
---|
135 | real(jprb), parameter :: tan_diffuse_angle_3d = Pi * 0.5_jprb |
---|
136 | |
---|
137 | ! The minimum cosine of solar zenith angle to consider for 3D |
---|
138 | ! effects, equivalent to one solar radius (0.2615 degrees) |
---|
139 | real(jprb), parameter :: min_mu0_3d = 0.004625_jprb |
---|
140 | |
---|
141 | ! Optical depth, single scattering albedo and asymmetry factor in |
---|
142 | ! each region and (except for asymmetry) at each g-point |
---|
143 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
144 | & :: od_region, ssa_region |
---|
145 | real(jprb), dimension(config%nregions) :: g_region |
---|
146 | |
---|
147 | ! Scattering optical depths of gases and clouds |
---|
148 | real(jprb) :: scat_od, scat_od_cloud |
---|
149 | |
---|
150 | ! The area fractions of each region |
---|
151 | real(jprb) :: region_fracs(1:config%nregions,nlev,istartcol:iendcol) |
---|
152 | |
---|
153 | ! The scaling used for the optical depth in the cloudy regions |
---|
154 | real(jprb) :: od_scaling(2:config%nregions,nlev,istartcol:iendcol) |
---|
155 | |
---|
156 | ! The length of the interface between regions jreg and jreg+1 per |
---|
157 | ! unit area of gridbox, equal to L_diff^ab in Hogan and Shonk |
---|
158 | ! (2013). When the first index is 3, this is the length of the |
---|
159 | ! interface between regions 3 and 1. This is actually the |
---|
160 | ! effective length oriented to a photon with random azimuth angle, |
---|
161 | ! so is the true edge length divided by pi. |
---|
162 | real(jprb) :: edge_length(3,nlev) |
---|
163 | |
---|
164 | ! Element i,j gives the rate of 3D transfer of diffuse/direct |
---|
165 | ! radiation from region i to region j, multiplied by the thickness |
---|
166 | ! of the layer in m |
---|
167 | real(jprb) :: transfer_rate_diffuse(config%nregions,config%nregions) |
---|
168 | real(jprb) :: transfer_rate_direct(config%nregions,config%nregions) |
---|
169 | |
---|
170 | ! Directional overlap matrices defined at all layer interfaces |
---|
171 | ! including top-of-atmosphere and the surface |
---|
172 | real(jprb), dimension(config%nregions,config%nregions,nlev+1, & |
---|
173 | & istartcol:iendcol) :: u_matrix, v_matrix |
---|
174 | |
---|
175 | ! Two-stream variables |
---|
176 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
177 | & :: gamma1, gamma2, gamma3 |
---|
178 | |
---|
179 | ! Matrix Gamma multiplied by the layer thickness z1, so units |
---|
180 | ! metres. After calling expm, this array contains the matrix |
---|
181 | ! exponential of the original. |
---|
182 | real(jprb) :: Gamma_z1(config%n_g_sw,3*config%nregions,3*config%nregions) |
---|
183 | |
---|
184 | ! Diffuse reflection and transmission matrices of each layer |
---|
185 | real(jprb), dimension(config%n_g_sw, config%nregions, & |
---|
186 | & config%nregions, nlev) :: reflectance, transmittance |
---|
187 | |
---|
188 | ! Clear-sky diffuse reflection and transmission matrices of each |
---|
189 | ! layer |
---|
190 | real(jprb), dimension(config%n_g_sw, nlev) :: ref_clear, trans_clear |
---|
191 | |
---|
192 | ! Matrices translating the direct flux entering the layer from |
---|
193 | ! above to the reflected radiation exiting upwards (ref_dir) and |
---|
194 | ! the scattered radiation exiting downwards (trans_dir_diff), |
---|
195 | ! along with the direct unscattered transmission matrix |
---|
196 | ! (trans_dir_dir). |
---|
197 | real(jprb), dimension(config%n_g_sw, config%nregions, config%nregions, nlev) & |
---|
198 | & :: ref_dir, trans_dir_diff, trans_dir_dir |
---|
199 | ! ...clear-sky equivalents |
---|
200 | real(jprb), dimension(config%n_g_sw, nlev) & |
---|
201 | & :: ref_dir_clear, trans_dir_diff_clear, trans_dir_dir_clear |
---|
202 | |
---|
203 | ! The fluxes downwelling from the bottom of the layer due to |
---|
204 | ! scattering by the direct beam within the layer |
---|
205 | real(jprb), dimension(config%n_g_sw, config%nregions) :: source_dn |
---|
206 | ! ...clear-sky equivalent |
---|
207 | real(jprb), dimension(config%n_g_sw) :: source_dn_clear |
---|
208 | |
---|
209 | ! The fluxes upwelling just above the base of a layer due to |
---|
210 | ! reflection of the direct downwelling beam; this is just used as |
---|
211 | ! a temporary variable |
---|
212 | real(jprb), dimension(config%n_g_sw, config%nregions) :: total_source |
---|
213 | |
---|
214 | ! Direct downwelling flux below and above an interface between |
---|
215 | ! layers into a plane perpendicular to the direction of the sun |
---|
216 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
217 | & :: direct_dn_below, direct_dn_above |
---|
218 | ! ...clear-sky equivalent (no distinction between "above/below") |
---|
219 | real(jprb), dimension(config%n_g_sw) :: direct_dn_clear |
---|
220 | |
---|
221 | ! Total albedo of the atmosphere/surface just above a layer |
---|
222 | ! interface with respect to downwelling diffuse and direct |
---|
223 | ! radiation at that interface, where level index = 1 corresponds |
---|
224 | ! to the top-of-atmosphere |
---|
225 | real(jprb), dimension(config%n_g_sw, config%nregions, & |
---|
226 | & config%nregions, nlev+1) :: total_albedo, total_albedo_direct |
---|
227 | ! ...clear-sky equivalent |
---|
228 | real(jprb), dimension(config%n_g_sw, nlev+1) & |
---|
229 | & :: total_albedo_clear, total_albedo_clear_direct |
---|
230 | |
---|
231 | ! As total_albedo, but just below a layer interface |
---|
232 | real(jprb), dimension(config%n_g_sw, config%nregions, config%nregions) & |
---|
233 | & :: total_albedo_below, total_albedo_below_direct |
---|
234 | |
---|
235 | ! Temporary array for applying adding method with entrapment to |
---|
236 | ! albedo matrices |
---|
237 | real(jprb), dimension(config%n_g_sw, config%nregions, config%nregions) & |
---|
238 | & :: albedo_part |
---|
239 | |
---|
240 | ! Horizontal migration distance (m) of reflected light |
---|
241 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
242 | & :: x_diffuse, x_direct |
---|
243 | ! Temporary variables when applying overlap rules |
---|
244 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
245 | & :: x_diffuse_above, x_direct_above |
---|
246 | |
---|
247 | real(jprb), dimension(config%n_g_sw, config%nregions, config%nregions) & |
---|
248 | & :: entrapment |
---|
249 | |
---|
250 | ! The following is used to store matrices of the form I-A*B that |
---|
251 | ! are used on the denominator of some expressions |
---|
252 | real(jprb) :: denominator(config%n_g_sw,config%nregions,config%nregions) |
---|
253 | |
---|
254 | ! Clear-sky equivalent, but actually its reciprocal to replace |
---|
255 | ! some divisions by multiplications |
---|
256 | real(jprb), dimension(config%n_g_sw) :: inv_denom_scalar |
---|
257 | |
---|
258 | ! Final step in working out how much transport between regions |
---|
259 | ! above occurs |
---|
260 | real(jprb), dimension(config%n_g_sw) :: fractal_factor |
---|
261 | |
---|
262 | ! Inverse of cloud effective size (m^-1) |
---|
263 | real(jprb) :: inv_effective_size |
---|
264 | |
---|
265 | ! Layer depth (m) |
---|
266 | real(jprb) :: dz, layer_depth(nlev) |
---|
267 | |
---|
268 | ! Upwelling and downwelling fluxes above/below layer interfaces |
---|
269 | real(jprb), dimension(config%n_g_sw, config%nregions) & |
---|
270 | & :: flux_up_above, flux_dn_above, flux_dn_below |
---|
271 | ! Clear-sky upwelling and downwelling fluxes (which don't |
---|
272 | ! distinguish between whether they are above/below a layer |
---|
273 | ! interface) |
---|
274 | real(jprb), dimension(config%n_g_sw) :: flux_up_clear, flux_dn_clear |
---|
275 | |
---|
276 | ! Index of top-most cloudy layer, or nlev+1 if no cloud |
---|
277 | integer :: i_cloud_top |
---|
278 | |
---|
279 | ! Keep a count of the number of calls to the two ways of computing |
---|
280 | ! reflectance/transmittance matrices |
---|
281 | integer :: n_calls_expm, n_calls_meador_weaver |
---|
282 | |
---|
283 | ! Identify clear-sky layers, with pseudo layers for outer space |
---|
284 | ! and below the ground, both treated as single-region clear skies |
---|
285 | logical :: is_clear_sky_layer(0:nlev+1) |
---|
286 | |
---|
287 | ! Used in computing rates of lateral radiation transfer |
---|
288 | real(jprb), parameter :: four_over_pi = 4.0_jprb / Pi |
---|
289 | |
---|
290 | ! Maximum entrapment coefficient |
---|
291 | real(jprb) :: max_entr |
---|
292 | |
---|
293 | real(jphook) :: hook_handle |
---|
294 | |
---|
295 | if (lhook) call dr_hook('radiation_spartacus_sw:solver_spartacus_sw',0,hook_handle) |
---|
296 | |
---|
297 | ! -------------------------------------------------------- |
---|
298 | ! Section 1: Prepare general variables and arrays |
---|
299 | ! -------------------------------------------------------- |
---|
300 | |
---|
301 | ! Copy array dimensions to local variables for convenience |
---|
302 | nreg = config%nregions |
---|
303 | ng = config%n_g_sw |
---|
304 | |
---|
305 | ! Reset count of number of calls to the two ways to compute |
---|
306 | ! reflectance/transmission matrices |
---|
307 | n_calls_expm = 0 |
---|
308 | n_calls_meador_weaver = 0 |
---|
309 | |
---|
310 | ! Compute the wavelength-independent region fractions and |
---|
311 | ! optical-depth scalings |
---|
312 | call calc_region_properties(nlev,nreg,istartcol,iendcol, & |
---|
313 | & config%i_cloud_pdf_shape == IPdfShapeGamma, & |
---|
314 | & cloud%fraction, cloud%fractional_std, region_fracs, & |
---|
315 | & od_scaling, config%cloud_fraction_threshold) |
---|
316 | |
---|
317 | ! Compute wavelength-independent overlap matrices u_matrix and v_matrix |
---|
318 | call calc_overlap_matrices(nlev,nreg,istartcol,iendcol, & |
---|
319 | & region_fracs, cloud%overlap_param, & |
---|
320 | & u_matrix, v_matrix, decorrelation_scaling=config%cloud_inhom_decorr_scaling, & |
---|
321 | & cloud_fraction_threshold=config%cloud_fraction_threshold, & |
---|
322 | & use_beta_overlap=config%use_beta_overlap, & |
---|
323 | & cloud_cover=flux%cloud_cover_sw) |
---|
324 | |
---|
325 | if (config%iverbose >= 3) then |
---|
326 | write(nulout,'(a)',advance='no') ' Processing columns' |
---|
327 | end if |
---|
328 | |
---|
329 | ! Main loop over columns |
---|
330 | do jcol = istartcol, iendcol |
---|
331 | ! -------------------------------------------------------- |
---|
332 | ! Section 2: Prepare column-specific variables and arrays |
---|
333 | ! -------------------------------------------------------- |
---|
334 | |
---|
335 | if (config%iverbose >= 3) then |
---|
336 | write(nulout,'(a)',advance='no') '.' |
---|
337 | end if |
---|
338 | |
---|
339 | ! Copy local cosine of the solar zenith angle |
---|
340 | mu0 = single_level%cos_sza(jcol) |
---|
341 | |
---|
342 | ! Skip profile if sun is too low in the sky |
---|
343 | if (mu0 < 1.0e-10_jprb) then |
---|
344 | flux%sw_dn(jcol,:) = 0.0_jprb |
---|
345 | flux%sw_up(jcol,:) = 0.0_jprb |
---|
346 | if (allocated(flux%sw_dn_direct)) then |
---|
347 | flux%sw_dn_direct(jcol,:) = 0.0_jprb |
---|
348 | end if |
---|
349 | if (config%do_clear) then |
---|
350 | flux%sw_dn_clear(jcol,:) = 0.0_jprb |
---|
351 | flux%sw_up_clear(jcol,:) = 0.0_jprb |
---|
352 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
353 | flux%sw_dn_direct_clear(jcol,:) = 0.0_jprb |
---|
354 | end if |
---|
355 | end if |
---|
356 | |
---|
357 | if (config%do_save_spectral_flux) then |
---|
358 | flux%sw_dn_band(:,jcol,:) = 0.0_jprb |
---|
359 | flux%sw_up_band(:,jcol,:) = 0.0_jprb |
---|
360 | if (allocated(flux%sw_dn_direct_band)) then |
---|
361 | flux%sw_dn_direct_band(:,jcol,:) = 0.0_jprb |
---|
362 | end if |
---|
363 | if (config%do_clear) then |
---|
364 | flux%sw_dn_clear_band(:,jcol,:) = 0.0_jprb |
---|
365 | flux%sw_up_clear_band(:,jcol,:) = 0.0_jprb |
---|
366 | if (allocated(flux%sw_dn_direct_clear_band)) then |
---|
367 | flux%sw_dn_direct_clear_band(:,jcol,:) = 0.0_jprb |
---|
368 | end if |
---|
369 | end if |
---|
370 | end if |
---|
371 | |
---|
372 | flux%sw_dn_diffuse_surf_g(:,jcol) = 0.0_jprb |
---|
373 | flux%sw_dn_direct_surf_g(:,jcol) = 0.0_jprb |
---|
374 | if (config%do_clear) then |
---|
375 | flux%sw_dn_diffuse_surf_clear_g(:,jcol) = 0.0_jprb |
---|
376 | flux%sw_dn_direct_surf_clear_g(:,jcol) = 0.0_jprb |
---|
377 | end if |
---|
378 | |
---|
379 | cycle |
---|
380 | end if ! sun is below the horizon |
---|
381 | |
---|
382 | ! At this point mu0 >= 1.0e-10 |
---|
383 | |
---|
384 | ! Used to compute rate of attenuation of direct solar beam |
---|
385 | ! through the atmosphere |
---|
386 | one_over_mu0 = 1.0_jprb / mu0 |
---|
387 | |
---|
388 | ! The rate at which direct radiation enters cloud sides is |
---|
389 | ! proportional to the tangent of the solar zenith angle, but |
---|
390 | ! this gets very large when the sun is low in the sky, in which |
---|
391 | ! case we limit it to one solar radius above the horizon |
---|
392 | if (mu0 < min_mu0_3d) then |
---|
393 | tan_sza = sqrt(1.0_jprb/(min_mu0_3d*min_mu0_3d) - 1.0_jprb) |
---|
394 | else if (one_over_mu0 > 1.0_jprb) then |
---|
395 | tan_sza = sqrt(one_over_mu0*one_over_mu0 - 1.0_jprb & |
---|
396 | & + config%overhead_sun_factor) |
---|
397 | else |
---|
398 | ! Just in case we get mu0 > 1... |
---|
399 | tan_sza = sqrt(config%overhead_sun_factor) |
---|
400 | end if |
---|
401 | |
---|
402 | ! Define which layers contain cloud; assume that |
---|
403 | ! cloud%crop_cloud_fraction has already been called |
---|
404 | is_clear_sky_layer = .true. |
---|
405 | i_cloud_top = nlev+1 |
---|
406 | do jlev = nlev,1,-1 |
---|
407 | if (cloud%fraction(jcol,jlev) > 0.0_jprb) then |
---|
408 | is_clear_sky_layer(jlev) = .false. |
---|
409 | i_cloud_top = jlev |
---|
410 | end if |
---|
411 | end do |
---|
412 | |
---|
413 | ! -------------------------------------------------------- |
---|
414 | ! Section 3: First loop over layers |
---|
415 | ! -------------------------------------------------------- |
---|
416 | ! In this section the reflectance, transmittance and sources |
---|
417 | ! are computed for each layer |
---|
418 | do jlev = 1,nlev ! Start at top-of-atmosphere |
---|
419 | ! -------------------------------------------------------- |
---|
420 | ! Section 3.1: Layer-specific initialization |
---|
421 | ! -------------------------------------------------------- |
---|
422 | |
---|
423 | ! Array-wise assignments |
---|
424 | gamma1 = 0.0_jprb |
---|
425 | gamma2 = 0.0_jprb |
---|
426 | gamma3 = 0.0_jprb |
---|
427 | Gamma_z1= 0.0_jprb |
---|
428 | transfer_rate_direct(:,:) = 0.0_jprb |
---|
429 | transfer_rate_diffuse(:,:) = 0.0_jprb |
---|
430 | edge_length(:,jlev) = 0.0_jprb |
---|
431 | |
---|
432 | ! The following is from the hydrostatic equation |
---|
433 | ! and ideal gas law: dz = dp * R * T / (p * g) |
---|
434 | layer_depth(jlev) = R_over_g & |
---|
435 | & * (thermodynamics%pressure_hl(jcol,jlev+1) & |
---|
436 | & - thermodynamics%pressure_hl(jcol,jlev)) & |
---|
437 | & * (thermodynamics%temperature_hl(jcol,jlev) & |
---|
438 | & + thermodynamics%temperature_hl(jcol,jlev+1)) & |
---|
439 | & / (thermodynamics%pressure_hl(jcol,jlev) & |
---|
440 | & + thermodynamics%pressure_hl(jcol,jlev+1)) |
---|
441 | |
---|
442 | ! -------------------------------------------------------- |
---|
443 | ! Section 3.2: Compute gamma variables |
---|
444 | ! -------------------------------------------------------- |
---|
445 | if (is_clear_sky_layer(jlev)) then |
---|
446 | ! --- Section 3.2a: Clear-sky case -------------------- |
---|
447 | |
---|
448 | nregactive = 1 ! Only region 1 (clear-sky) is active |
---|
449 | |
---|
450 | ! Copy optical depth and single-scattering albedo of |
---|
451 | ! clear-sky region |
---|
452 | od_region(1:ng,1) = od(1:ng,jlev,jcol) |
---|
453 | ssa_region(1:ng,1) = ssa(1:ng,jlev,jcol) |
---|
454 | |
---|
455 | ! Compute two-stream variables gamma1-gamma3 (gamma4=1-gamma3) |
---|
456 | call calc_two_stream_gammas_sw(ng, & |
---|
457 | & mu0, ssa(1:ng,jlev,jcol), g(1:ng,jlev,jcol), & |
---|
458 | & gamma1(1:ng,1), gamma2(1:ng,1), gamma3(1:ng,1)) |
---|
459 | |
---|
460 | if (config%use_expm_everywhere) then |
---|
461 | ! Treat cloud-free layer as 3D: the matrix exponential |
---|
462 | ! "expm" is used as much as possible (provided the |
---|
463 | ! optical depth is not too high). ng3D is initially |
---|
464 | ! set to the total number of g-points, then the |
---|
465 | ! clear-sky optical depths are searched until a value |
---|
466 | ! exceeds the threshold for treating 3D effects, and |
---|
467 | ! if found it is reset to that. Note that we are |
---|
468 | ! assuming that the g-points have been reordered in |
---|
469 | ! approximate order of gas optical depth. |
---|
470 | ng3D = ng |
---|
471 | do jg = 1, ng |
---|
472 | if (od_region(jg,1) > config%max_gas_od_3D) then |
---|
473 | ng3D = jg-1 |
---|
474 | exit |
---|
475 | end if |
---|
476 | end do |
---|
477 | else |
---|
478 | ! Otherwise treat cloud-free layer using the classical |
---|
479 | ! Meador & Weaver formulae for all g-points |
---|
480 | ng3D = 0 |
---|
481 | end if |
---|
482 | |
---|
483 | |
---|
484 | else |
---|
485 | ! --- Section 3.2b: Cloudy case ----------------------- |
---|
486 | |
---|
487 | ! Default number of g-points to treat with |
---|
488 | ! matrix-exponential scheme |
---|
489 | if (config%use_expm_everywhere) then |
---|
490 | ng3D = ng ! All g-points |
---|
491 | else |
---|
492 | ng3D = 0 ! No g-points |
---|
493 | end if |
---|
494 | |
---|
495 | if (config%do_3d_effects .and. & |
---|
496 | & allocated(cloud%inv_cloud_effective_size) .and. & |
---|
497 | & .not. (nreg == 2 .and. cloud%fraction(jcol,jlev) & |
---|
498 | & > 1.0-config%cloud_fraction_threshold)) then |
---|
499 | if (cloud%inv_cloud_effective_size(jcol,jlev) > 0.0_jprb) then |
---|
500 | ! 3D effects are only simulated if |
---|
501 | ! inv_cloud_effective_size is defined and greater |
---|
502 | ! than zero |
---|
503 | ng3D = ng |
---|
504 | |
---|
505 | ! Depth of current layer |
---|
506 | dz = layer_depth(jlev) |
---|
507 | |
---|
508 | ! Compute cloud edge length per unit area of gridbox |
---|
509 | ! from rearranging Hogan & Shonk (2013) Eq. 45, but |
---|
510 | ! adding a factor of (1-frac) so that a region that |
---|
511 | ! fully occupies the gridbox (frac=1) has an edge of |
---|
512 | ! zero. We can therefore use the fraction of clear sky, |
---|
513 | ! region_fracs(1,jlev,jcol) for convenience instead. The |
---|
514 | ! pi on the denominator means that this is actually edge |
---|
515 | ! length with respect to a light ray with a random |
---|
516 | ! azimuthal direction. |
---|
517 | edge_length(1,jlev) = four_over_pi & |
---|
518 | & * region_fracs(1,jlev,jcol)*(1.0_jprb-region_fracs(1,jlev,jcol)) & |
---|
519 | & * min(cloud%inv_cloud_effective_size(jcol,jlev), & |
---|
520 | & 1.0_jprb / config%min_cloud_effective_size) |
---|
521 | if (nreg > 2) then |
---|
522 | ! The corresponding edge length between the two cloudy |
---|
523 | ! regions is computed in the same way but treating the |
---|
524 | ! optically denser of the two regions (region 3) as |
---|
525 | ! the cloud; note that the fraction of this region may |
---|
526 | ! be less than that of the optically less dense of the |
---|
527 | ! two regions (region 2). For increased flexibility, |
---|
528 | ! the user may specify the effective size of |
---|
529 | ! inhomogeneities separately from the cloud effective |
---|
530 | ! size. |
---|
531 | if (allocated(cloud%inv_inhom_effective_size)) then |
---|
532 | edge_length(2,jlev) = four_over_pi & |
---|
533 | & * region_fracs(3,jlev,jcol)*(1.0_jprb-region_fracs(3,jlev,jcol)) & |
---|
534 | & * min(cloud%inv_inhom_effective_size(jcol,jlev), & |
---|
535 | & 1.0_jprb / config%min_cloud_effective_size) |
---|
536 | else |
---|
537 | edge_length(2,jlev) = four_over_pi & |
---|
538 | & * region_fracs(3,jlev,jcol)*(1.0_jprb-region_fracs(3,jlev,jcol)) & |
---|
539 | & * min(cloud%inv_cloud_effective_size(jcol,jlev), & |
---|
540 | & 1.0_jprb / config%min_cloud_effective_size) |
---|
541 | end if |
---|
542 | |
---|
543 | ! In the case of three regions, some of the cloud |
---|
544 | ! boundary may go directly to the "thick" region |
---|
545 | if (config%clear_to_thick_fraction > 0.0_jprb) then |
---|
546 | edge_length(3,jlev) = config%clear_to_thick_fraction & |
---|
547 | & * min(edge_length(1,jlev), edge_length(2,jlev)) |
---|
548 | edge_length(1,jlev) = edge_length(1,jlev) - edge_length(3,jlev) |
---|
549 | edge_length(2,jlev) = edge_length(2,jlev) - edge_length(3,jlev) |
---|
550 | else |
---|
551 | edge_length(3,jlev) = 0.0_jprb |
---|
552 | end if |
---|
553 | end if |
---|
554 | |
---|
555 | do jreg = 1, nreg-1 |
---|
556 | ! Compute lateral transfer rates from region jreg to |
---|
557 | ! jreg+1 following Hogan & Shonk (2013) Eq. 47, but |
---|
558 | ! multiplied by dz because the transfer rate is |
---|
559 | ! vertically integrated across the depth of the layer |
---|
560 | if (region_fracs(jreg,jlev,jcol) > epsilon(1.0_jprb)) then |
---|
561 | transfer_rate_direct(jreg,jreg+1) = dz & |
---|
562 | & * edge_length(jreg,jlev) * tan_sza / region_fracs(jreg,jlev,jcol) |
---|
563 | transfer_rate_diffuse(jreg,jreg+1) = dz & |
---|
564 | & * edge_length(jreg,jlev) & |
---|
565 | & * tan_diffuse_angle_3d / region_fracs(jreg,jlev,jcol) |
---|
566 | end if |
---|
567 | ! Compute transfer rates from region jreg+1 to |
---|
568 | ! jreg |
---|
569 | if (region_fracs(jreg+1,jlev,jcol) > epsilon(1.0_jprb)) then |
---|
570 | transfer_rate_direct(jreg+1,jreg) = dz & |
---|
571 | & * edge_length(jreg,jlev) & |
---|
572 | & * tan_sza / region_fracs(jreg+1,jlev,jcol) |
---|
573 | transfer_rate_diffuse(jreg+1,jreg) = dz & |
---|
574 | & * edge_length(jreg,jlev) & |
---|
575 | & * tan_diffuse_angle_3d / region_fracs(jreg+1,jlev,jcol) |
---|
576 | end if |
---|
577 | end do |
---|
578 | |
---|
579 | ! Compute transfer rates directly between regions 1 and |
---|
580 | ! 3 |
---|
581 | if (edge_length(3,jlev) > 0.0_jprb) then |
---|
582 | if (region_fracs(1,jlev,jcol) > epsilon(1.0_jprb)) then |
---|
583 | transfer_rate_direct(1,3) = dz & |
---|
584 | & * edge_length(3,jlev) * tan_sza / region_fracs(1,jlev,jcol) |
---|
585 | transfer_rate_diffuse(1,3) = dz & |
---|
586 | & * edge_length(3,jlev) & |
---|
587 | & * tan_diffuse_angle_3d / region_fracs(1,jlev,jcol) |
---|
588 | end if |
---|
589 | if (region_fracs(3,jlev,jcol) > epsilon(1.0_jprb)) then |
---|
590 | transfer_rate_direct(3,1) = dz & |
---|
591 | & * edge_length(3,jlev) * tan_sza / region_fracs(3,jlev,jcol) |
---|
592 | transfer_rate_diffuse(3,1) = dz & |
---|
593 | & * edge_length(3,jlev) & |
---|
594 | & * tan_diffuse_angle_3d / region_fracs(3,jlev,jcol) |
---|
595 | end if |
---|
596 | end if |
---|
597 | |
---|
598 | ! Don't allow the transfer rate out of a region to be |
---|
599 | ! equivalent to a loss of exp(-10) through the layer |
---|
600 | where (transfer_rate_direct > config%max_3d_transfer_rate) |
---|
601 | transfer_rate_direct = config%max_3d_transfer_rate |
---|
602 | end where |
---|
603 | where (transfer_rate_diffuse > config%max_3d_transfer_rate) |
---|
604 | transfer_rate_diffuse = config%max_3d_transfer_rate |
---|
605 | end where |
---|
606 | |
---|
607 | end if ! Cloud has edge length required for 3D effects |
---|
608 | end if ! Include 3D effects |
---|
609 | |
---|
610 | ! In a cloudy layer the number of active regions equals |
---|
611 | ! the number of regions |
---|
612 | nregactive = nreg |
---|
613 | |
---|
614 | ! Compute scattering properties of the regions at each |
---|
615 | ! g-point, mapping from the cloud properties |
---|
616 | ! defined in each band. |
---|
617 | do jg = 1,ng |
---|
618 | ! Mapping from g-point to band |
---|
619 | iband = config%i_band_from_reordered_g_sw(jg) |
---|
620 | |
---|
621 | ! Scattering optical depth of clear-sky region |
---|
622 | scat_od = od(jg,jlev,jcol)*ssa(jg,jlev,jcol) |
---|
623 | |
---|
624 | ! Scattering properties of clear-sky regions copied |
---|
625 | ! over |
---|
626 | od_region(jg,1) = od(jg, jlev, jcol) |
---|
627 | ssa_region(jg,1) = ssa(jg, jlev, jcol) |
---|
628 | g_region(1) = g(jg, jlev, jcol) |
---|
629 | |
---|
630 | ! Loop over cloudy regions |
---|
631 | do jreg = 2,nreg |
---|
632 | scat_od_cloud = od_cloud(iband,jlev,jcol) & |
---|
633 | & * ssa_cloud(iband,jlev,jcol)*od_scaling(jreg,jlev,jcol) |
---|
634 | ! Add scaled cloud optical depth to clear-sky value |
---|
635 | od_region(jg,jreg) = od(jg,jlev,jcol) & |
---|
636 | & + od_cloud(iband,jlev,jcol)*od_scaling(jreg,jlev,jcol) |
---|
637 | ! Compute single-scattering albedo and asymmetry |
---|
638 | ! factor of gas-cloud combination |
---|
639 | ssa_region(jg,jreg) = (scat_od+scat_od_cloud) & |
---|
640 | & / od_region(jg,jreg) |
---|
641 | g_region(jreg) = (scat_od*g(jg,jlev,jcol) & |
---|
642 | & + scat_od_cloud * g_cloud(iband,jlev,jcol)) & |
---|
643 | & / (scat_od + scat_od_cloud) |
---|
644 | |
---|
645 | ! Apply maximum cloud optical depth for stability in the |
---|
646 | ! 3D case |
---|
647 | if (od_region(jg,jreg) > config%max_cloud_od) then |
---|
648 | od_region(jg,jreg) = config%max_cloud_od |
---|
649 | end if |
---|
650 | |
---|
651 | end do |
---|
652 | |
---|
653 | ! Calculate two-stream variables gamma1-gamma3 of all |
---|
654 | ! regions at once |
---|
655 | call calc_two_stream_gammas_sw(nreg, & |
---|
656 | & mu0, ssa_region(jg,:), g_region, & |
---|
657 | & gamma1(jg,:), gamma2(jg,:), gamma3(jg,:)) |
---|
658 | |
---|
659 | ! Loop is in order of g-points with typically |
---|
660 | ! increasing optical depth: if optical depth of |
---|
661 | ! clear-sky region exceeds a threshold then turn off |
---|
662 | ! 3D effects for any further g-points |
---|
663 | if (ng3D == ng & |
---|
664 | & .and. od_region(jg,1) > config%max_gas_od_3D) then |
---|
665 | ng3D = jg-1 |
---|
666 | end if |
---|
667 | end do ! Loop over g points |
---|
668 | end if ! Cloudy level |
---|
669 | |
---|
670 | ! -------------------------------------------------------- |
---|
671 | ! Section 3.3: Compute reflection, transmission and emission |
---|
672 | ! -------------------------------------------------------- |
---|
673 | if (ng3D > 0) then |
---|
674 | ! --- Section 3.3a: g-points with 3D effects ---------- |
---|
675 | |
---|
676 | ! 3D effects need to be represented in "ng3D" of the g |
---|
677 | ! points. This is done by creating ng3D square matrices |
---|
678 | ! each of dimension 3*nreg by 3*nreg, computing the matrix |
---|
679 | ! exponential, then computing the various |
---|
680 | ! transmission/reflectance matrices from that. |
---|
681 | do jreg = 1,nregactive |
---|
682 | ! Write the diagonal elements of -Gamma1*z1 |
---|
683 | Gamma_z1(1:ng3D,jreg,jreg) & |
---|
684 | & = od_region(1:ng3D,jreg)*gamma1(1:ng3D,jreg) |
---|
685 | ! Write the diagonal elements of +Gamma2*z1 |
---|
686 | Gamma_z1(1:ng3D,jreg+nreg,jreg) & |
---|
687 | & = od_region(1:ng3D,jreg)*gamma2(1:ng3D,jreg) |
---|
688 | ! Write the diagonal elements of -Gamma3*z1 |
---|
689 | Gamma_z1(1:ng3D,jreg,jreg+2*nreg) & |
---|
690 | & = -od_region(1:ng3D,jreg)*ssa_region(1:ng3D,jreg) & |
---|
691 | & * gamma3(1:ng3D,jreg) |
---|
692 | |
---|
693 | ! Write the diagonal elements of +Gamma4*z1 |
---|
694 | Gamma_z1(1:ng3D,jreg+nreg,jreg+2*nreg) & |
---|
695 | & = od_region(1:ng3D,jreg)*ssa_region(1:ng3D,jreg) & |
---|
696 | & * (1.0_jprb - gamma3(1:ng3D,jreg)) |
---|
697 | |
---|
698 | ! Write the diagonal elements of +Gamma0*z1 |
---|
699 | Gamma_z1(1:ng3D,jreg+2*nreg,jreg+2*nreg) & |
---|
700 | & = -od_region(1:ng3D,jreg)*one_over_mu0 |
---|
701 | end do |
---|
702 | |
---|
703 | do jreg = 1,nregactive-1 |
---|
704 | ! Write the elements of -Gamma1*z1 concerned with 3D |
---|
705 | ! transport |
---|
706 | Gamma_z1(1:ng3D,jreg,jreg) = Gamma_z1(1:ng3D,jreg,jreg) & |
---|
707 | & + transfer_rate_diffuse(jreg,jreg+1) |
---|
708 | Gamma_z1(1:ng3D,jreg+1,jreg+1) = Gamma_z1(1:ng3D,jreg+1,jreg+1) & |
---|
709 | & + transfer_rate_diffuse(jreg+1,jreg) |
---|
710 | Gamma_z1(1:ng3D,jreg+1,jreg) = -transfer_rate_diffuse(jreg,jreg+1) |
---|
711 | Gamma_z1(1:ng3D,jreg,jreg+1) = -transfer_rate_diffuse(jreg+1,jreg) |
---|
712 | ! Write the elements of +Gamma0*z1 concerned with 3D |
---|
713 | ! transport |
---|
714 | Gamma_z1(1:ng3D,jreg+2*nreg,jreg+2*nreg) & |
---|
715 | & = Gamma_z1(1:ng3D,jreg+2*nreg,jreg+2*nreg) & |
---|
716 | & - transfer_rate_direct(jreg,jreg+1) |
---|
717 | Gamma_z1(1:ng3D,jreg+2*nreg+1,jreg+2*nreg+1) & |
---|
718 | & = Gamma_z1(1:ng3D,jreg+2*nreg+1,jreg+2*nreg+1) & |
---|
719 | & - transfer_rate_direct(jreg+1,jreg) |
---|
720 | Gamma_z1(1:ng3D,jreg+2*nreg+1,jreg+2*nreg) & |
---|
721 | & = transfer_rate_direct(jreg,jreg+1) |
---|
722 | Gamma_z1(1:ng3D,jreg+2*nreg,jreg+2*nreg+1) & |
---|
723 | & = transfer_rate_direct(jreg+1,jreg) |
---|
724 | end do |
---|
725 | |
---|
726 | ! Possible flow between regions a and c |
---|
727 | if (edge_length(3,jlev) > 0.0_jprb) then |
---|
728 | ! Diffuse transport |
---|
729 | Gamma_z1(1:ng3D,1,1) = Gamma_z1(1:ng3D,1,1) & |
---|
730 | & + transfer_rate_diffuse(1,3) |
---|
731 | Gamma_z1(1:ng3D,3,3) = Gamma_z1(1:ng3D,3,3) & |
---|
732 | & + transfer_rate_diffuse(3,1) |
---|
733 | Gamma_z1(1:ng3D,3,1) = -transfer_rate_diffuse(1,3) |
---|
734 | Gamma_z1(1:ng3D,1,3) = -transfer_rate_diffuse(3,1) |
---|
735 | ! Direct transport |
---|
736 | Gamma_z1(1:ng3D,1+2*nreg,1+2*nreg) = Gamma_z1(1:ng3D,1+2*nreg,1+2*nreg) & |
---|
737 | & - transfer_rate_direct(1,3) |
---|
738 | Gamma_z1(1:ng3D,3+2*nreg,3+2*nreg) = Gamma_z1(1:ng3D,3+2*nreg,3+2*nreg) & |
---|
739 | & - transfer_rate_direct(3,1) |
---|
740 | Gamma_z1(1:ng3D,3+2*nreg,1+2*nreg) = transfer_rate_direct(1,3) |
---|
741 | Gamma_z1(1:ng3D,1+2*nreg,3+2*nreg) = transfer_rate_direct(3,1) |
---|
742 | end if |
---|
743 | |
---|
744 | ! Copy Gamma1*z1 |
---|
745 | Gamma_z1(1:ng3D,nreg+1:nreg+nregactive,nreg+1:nreg+nregactive) & |
---|
746 | & = -Gamma_z1(1:ng3D,1:nregactive,1:nregactive) |
---|
747 | ! Copy Gamma2*z1 |
---|
748 | Gamma_z1(1:ng3D,1:nregactive,nreg+1:nreg+nregactive) & |
---|
749 | & = -Gamma_z1(1:ng3D,nreg+1:nreg+nregactive,1:nregactive) |
---|
750 | |
---|
751 | ! Compute the matrix exponential of Gamma_z1, returning the |
---|
752 | ! result in-place |
---|
753 | call expm(ng, ng3D, 3*nreg, Gamma_z1, IMatrixPatternShortwave) |
---|
754 | |
---|
755 | ! Update count of expm calls |
---|
756 | n_calls_expm = n_calls_expm + ng3D |
---|
757 | |
---|
758 | ! Direct transmission matrix |
---|
759 | trans_dir_dir(1:ng3D,:,:,jlev) = min(1.0_jprb,max(0.0_jprb, & |
---|
760 | & Gamma_z1(1:ng3D,2*nreg+1:3*nreg, 2*nreg+1:3*nreg))) |
---|
761 | ! Diffuse reflectance matrix; security on negative values |
---|
762 | ! necessary occasionally for very low cloud fraction and very high |
---|
763 | ! in-cloud optical depth |
---|
764 | reflectance(1:ng3D,:,:,jlev) = min(1.0_jprb,max(0.0_jprb, & |
---|
765 | & -solve_mat(ng,ng3D,nreg,Gamma_z1(1:ng3D,1:nreg,1:nreg), & |
---|
766 | & Gamma_z1(1:ng3D,1:nreg,nreg+1:2*nreg)))) |
---|
767 | ! Diffuse transmission matrix |
---|
768 | transmittance(1:ng3D,:,:,jlev) = min(1.0_jprb,max(0.0_jprb, & |
---|
769 | & mat_x_mat(ng,ng3D,nreg,Gamma_z1(1:ng3D,nreg+1:2*nreg,1:nreg), & |
---|
770 | & reflectance(1:ng3D,:,:,jlev)) & |
---|
771 | & + Gamma_z1(1:ng3D,nreg+1:2*nreg,nreg+1:2*nreg))) |
---|
772 | ! Transfer matrix between downward direct and upward |
---|
773 | ! diffuse |
---|
774 | ref_dir(1:ng3D,:,:,jlev) = min(mu0,max(0.0_jprb, & |
---|
775 | & -solve_mat(ng,ng3D,nreg,Gamma_z1(1:ng3D,1:nreg,1:nreg), & |
---|
776 | & Gamma_z1(1:ng3D,1:nreg,2*nreg+1:3*nreg)))) |
---|
777 | ! Transfer matrix between downward direct and downward |
---|
778 | ! diffuse in layer interface below. Include correction for |
---|
779 | ! trans_dir_diff out of plausible bounds (note that Meador & |
---|
780 | ! Weaver has the same correction in radiation_two_stream.F90 |
---|
781 | ! - this is not just an expm thing) |
---|
782 | trans_dir_diff(1:ng3D,:,:,jlev) = min(mu0,max(0.0_jprb, & |
---|
783 | & mat_x_mat(ng,ng3D,nreg,Gamma_z1(1:ng3D,nreg+1:2*nreg,1:nreg), & |
---|
784 | & ref_dir(1:ng3D,:,:,jlev)) & |
---|
785 | & + Gamma_z1(1:ng3D,nreg+1:2*nreg,2*nreg+1:3*nreg))) |
---|
786 | |
---|
787 | end if ! we are treating 3D effects for some g points |
---|
788 | |
---|
789 | ! --- Section 3.3b: g-points without 3D effects ---------- |
---|
790 | |
---|
791 | ! Compute reflectance, transmittance and associated terms for |
---|
792 | ! clear skies, using the Meador-Weaver formulas |
---|
793 | call calc_reflectance_transmittance_sw(ng, & |
---|
794 | & mu0, od_region(1:ng,1), ssa_region(1:ng,1), & |
---|
795 | & gamma1(1:ng,1), gamma2(1:ng,1), gamma3(1:ng,1), & |
---|
796 | & ref_clear(1:ng,jlev), trans_clear(1:ng,jlev), & |
---|
797 | & ref_dir_clear(1:ng,jlev), trans_dir_diff_clear(1:ng,jlev), & |
---|
798 | & trans_dir_dir_clear(1:ng,jlev) ) |
---|
799 | |
---|
800 | n_calls_meador_weaver = n_calls_meador_weaver + ng |
---|
801 | |
---|
802 | if (ng3D < ng) then |
---|
803 | ! Some of the g points are to be treated using the |
---|
804 | ! conventional plane-parallel method. First zero the |
---|
805 | ! relevant parts of the matrices |
---|
806 | trans_dir_dir (ng3D+1:ng,:,:,jlev) = 0.0_jprb |
---|
807 | reflectance (ng3D+1:ng,:,:,jlev) = 0.0_jprb |
---|
808 | transmittance (ng3D+1:ng,:,:,jlev) = 0.0_jprb |
---|
809 | ref_dir (ng3D+1:ng,:,:,jlev) = 0.0_jprb |
---|
810 | trans_dir_diff(ng3D+1:ng,:,:,jlev) = 0.0_jprb |
---|
811 | |
---|
812 | ! Since there is no lateral transport, the clear-sky parts |
---|
813 | ! of the arrays can be copied from the clear-sky arrays |
---|
814 | trans_dir_dir (ng3D+1:ng,1,1,jlev) = trans_dir_dir_clear (ng3D+1:ng,jlev) |
---|
815 | reflectance (ng3D+1:ng,1,1,jlev) = ref_clear (ng3D+1:ng,jlev) |
---|
816 | transmittance (ng3D+1:ng,1,1,jlev) = trans_clear (ng3D+1:ng,jlev) |
---|
817 | ref_dir (ng3D+1:ng,1,1,jlev) = ref_dir_clear (ng3D+1:ng,jlev) |
---|
818 | trans_dir_diff(ng3D+1:ng,1,1,jlev) = trans_dir_diff_clear(ng3D+1:ng,jlev) |
---|
819 | |
---|
820 | ! Compute reflectance, transmittance and associated terms |
---|
821 | ! for each cloudy region, using the Meador-Weaver formulas |
---|
822 | do jreg = 2, nregactive |
---|
823 | call calc_reflectance_transmittance_sw(ng-ng3D, & |
---|
824 | & mu0, & |
---|
825 | & od_region(ng3D+1:ng,jreg), ssa_region(ng3D+1:ng,jreg), & |
---|
826 | & gamma1(ng3D+1:ng,jreg), gamma2(ng3D+1:ng,jreg), & |
---|
827 | & gamma3(ng3D+1:ng,jreg), & |
---|
828 | & reflectance(ng3D+1:ng,jreg,jreg,jlev), & |
---|
829 | & transmittance(ng3D+1:ng,jreg,jreg,jlev), & |
---|
830 | & ref_dir(ng3D+1:ng,jreg,jreg,jlev), & |
---|
831 | & trans_dir_diff(ng3D+1:ng,jreg,jreg,jlev), & |
---|
832 | & trans_dir_dir(ng3D+1:ng,jreg,jreg,jlev) ) |
---|
833 | end do |
---|
834 | n_calls_meador_weaver & |
---|
835 | & = n_calls_meador_weaver + (ng-ng3D)*(nregactive-1) |
---|
836 | end if |
---|
837 | |
---|
838 | end do ! Loop over levels |
---|
839 | |
---|
840 | ! -------------------------------------------------------- |
---|
841 | ! Section 4: Compute total albedos |
---|
842 | ! -------------------------------------------------------- |
---|
843 | |
---|
844 | total_albedo(:,:,:,:) = 0.0_jprb |
---|
845 | total_albedo_direct(:,:,:,:) = 0.0_jprb |
---|
846 | |
---|
847 | if (config%do_clear) then |
---|
848 | total_albedo_clear(:,:) = 0.0_jprb |
---|
849 | total_albedo_clear_direct(:,:) = 0.0_jprb |
---|
850 | end if |
---|
851 | |
---|
852 | ! Calculate the upwelling radiation scattered from the direct |
---|
853 | ! beam incident on the surface, and copy the surface albedo |
---|
854 | ! into total_albedo |
---|
855 | do jreg = 1,nreg |
---|
856 | do jg = 1,ng |
---|
857 | total_albedo(jg,jreg,jreg,nlev+1) = albedo_diffuse(jg,jcol) |
---|
858 | total_albedo_direct(jg,jreg,jreg,nlev+1) & |
---|
859 | & = mu0 * albedo_direct(jg,jcol) |
---|
860 | end do |
---|
861 | end do |
---|
862 | |
---|
863 | if (config%do_clear) then |
---|
864 | ! Surface albedo is the same |
---|
865 | total_albedo_clear(1:ng,nlev+1) = total_albedo(1:ng,1,1,nlev+1) |
---|
866 | total_albedo_clear_direct(1:ng,nlev+1) & |
---|
867 | & = total_albedo_direct(1:ng,1,1,nlev+1) |
---|
868 | end if |
---|
869 | |
---|
870 | ! Horizontal migration distances of reflected radiation at the |
---|
871 | ! surface are zero |
---|
872 | x_diffuse = 0.0_jprb |
---|
873 | x_direct = 0.0_jprb |
---|
874 | |
---|
875 | ! Work back up through the atmosphere computing the total albedo |
---|
876 | ! of the atmosphere below that point using the adding method |
---|
877 | do jlev = nlev,1,-1 |
---|
878 | |
---|
879 | ! -------------------------------------------------------- |
---|
880 | ! Section 4.1: Adding method |
---|
881 | ! -------------------------------------------------------- |
---|
882 | |
---|
883 | if (config%do_clear) then |
---|
884 | ! Use adding method for clear-sky arrays; note that there |
---|
885 | ! is no need to consider "above" and "below" quantities |
---|
886 | ! since with no cloud overlap to worry about, these are |
---|
887 | ! the same |
---|
888 | inv_denom_scalar(:) = 1.0_jprb & |
---|
889 | & / (1.0_jprb - total_albedo_clear(:,jlev+1)*ref_clear(:,jlev)) |
---|
890 | total_albedo_clear(:,jlev) = ref_clear(:,jlev) & |
---|
891 | & + trans_clear(:,jlev)*trans_clear(:,jlev)*total_albedo_clear(:,jlev+1) & |
---|
892 | & * inv_denom_scalar(:) |
---|
893 | total_albedo_clear_direct(:,jlev) = ref_dir_clear(:,jlev) & |
---|
894 | & + (trans_dir_dir_clear(:,jlev) * total_albedo_clear_direct(:,jlev+1) & |
---|
895 | & +trans_dir_diff_clear(:,jlev) * total_albedo_clear(:,jlev+1)) & |
---|
896 | & * trans_clear(:,jlev) * inv_denom_scalar(:) |
---|
897 | end if |
---|
898 | |
---|
899 | if (is_clear_sky_layer(jlev)) then |
---|
900 | ! Clear-sky layer: use scalar adding method |
---|
901 | inv_denom_scalar(:) = 1.0_jprb & |
---|
902 | & / (1.0_jprb - total_albedo(:,1,1,jlev+1)*reflectance(:,1,1,jlev)) |
---|
903 | total_albedo_below = 0.0_jprb |
---|
904 | total_albedo_below(:,1,1) = reflectance(:,1,1,jlev) & |
---|
905 | & + transmittance(:,1,1,jlev) * transmittance(:,1,1,jlev) & |
---|
906 | & * total_albedo(:,1,1,jlev+1) * inv_denom_scalar(:) |
---|
907 | total_albedo_below_direct = 0.0_jprb |
---|
908 | total_albedo_below_direct(:,1,1) = ref_dir(:,1,1,jlev) & |
---|
909 | & + (trans_dir_dir(:,1,1,jlev)*total_albedo_direct(:,1,1,jlev+1) & |
---|
910 | & +trans_dir_diff(:,1,1,jlev)*total_albedo(:,1,1,jlev+1)) & |
---|
911 | & * transmittance(:,1,1,jlev) * inv_denom_scalar(:) |
---|
912 | else |
---|
913 | ! Cloudy layer: use matrix adding method |
---|
914 | denominator = identity_minus_mat_x_mat(ng,ng,nreg, & |
---|
915 | & total_albedo(:,:,:,jlev+1), reflectance(:,:,:,jlev)) |
---|
916 | total_albedo_below = reflectance(:,:,:,jlev) & |
---|
917 | & + mat_x_mat(ng,ng,nreg,transmittance(:,:,:,jlev), & |
---|
918 | & solve_mat(ng,ng,nreg,denominator, & |
---|
919 | & mat_x_mat(ng,ng,nreg,total_albedo(:,:,:,jlev+1), & |
---|
920 | & transmittance(:,:,:,jlev)))) |
---|
921 | total_albedo_below_direct = ref_dir(:,:,:,jlev) & |
---|
922 | & + mat_x_mat(ng,ng,nreg,transmittance(:,:,:,jlev), & |
---|
923 | & solve_mat(ng,ng,nreg,denominator, & |
---|
924 | & mat_x_mat(ng,ng,nreg,total_albedo_direct(:,:,:,jlev+1), & |
---|
925 | & trans_dir_dir(:,:,:,jlev)) & |
---|
926 | & +mat_x_mat(ng,ng,nreg,total_albedo(:,:,:,jlev+1), & |
---|
927 | & trans_dir_diff(:,:,:,jlev)))) |
---|
928 | end if |
---|
929 | |
---|
930 | ! -------------------------------------------------------- |
---|
931 | ! Section 4.2: Overlap and entrapment |
---|
932 | ! -------------------------------------------------------- |
---|
933 | |
---|
934 | #ifndef PRINT_ENTRAPMENT_DATA |
---|
935 | if ((config%i_3d_sw_entrapment == IEntrapmentExplicitNonFractal & |
---|
936 | & .or. config%i_3d_sw_entrapment == IEntrapmentExplicit) & |
---|
937 | & .and. jlev >= i_cloud_top) then |
---|
938 | #else |
---|
939 | if (config%i_3d_sw_entrapment == IEntrapmentExplicitNonFractal & |
---|
940 | & .or. config%i_3d_sw_entrapment == IEntrapmentExplicit) then |
---|
941 | #endif |
---|
942 | ! "Explicit entrapment": we have the horizontal migration |
---|
943 | ! distances just above the base of the layer, and need to |
---|
944 | ! step them to just below the top of the same layer |
---|
945 | call step_migrations(ng, nreg, cloud%fraction(jcol,jlev), & |
---|
946 | & layer_depth(jlev), tan_diffuse_angle_3d, tan_sza, & |
---|
947 | & reflectance(:,:,:,jlev), transmittance(:,:,:,jlev), & |
---|
948 | & ref_dir(:,:,:,jlev), trans_dir_dir(:,:,:,jlev), & |
---|
949 | & trans_dir_diff(:,:,:,jlev), total_albedo(:,:,:,jlev+1), & |
---|
950 | & total_albedo_direct(:,:,:,jlev+1), & |
---|
951 | & x_diffuse, x_direct) |
---|
952 | |
---|
953 | #ifdef PRINT_ENTRAPMENT_DATA |
---|
954 | ! Write out for later analysis: these are the entrapment |
---|
955 | ! statistics at the top of layer "jlev" |
---|
956 | ! Note that number of scattering events is now not computed, |
---|
957 | ! so print "1.0" |
---|
958 | if (nreg == 2) then |
---|
959 | write(101,'(i4,i4,6e14.6)') jcol, jlev, & |
---|
960 | & x_direct(1,:), x_diffuse(1,:), x_direct(1,:)*0.0_jprb+1.0_jprb |
---|
961 | else |
---|
962 | write(101,'(i4,i4,9e14.6)') jcol, jlev, & |
---|
963 | & x_direct(1,1:3), x_diffuse(1,1:3), 1.0_jprb,1.0_jprb,1.0_jprb |
---|
964 | end if |
---|
965 | #endif |
---|
966 | |
---|
967 | end if |
---|
968 | |
---|
969 | ! Account for cloud overlap when converting albedo and source |
---|
970 | ! below a layer interface to the equivalent values just above |
---|
971 | if (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev-1)) then |
---|
972 | ! If both layers are cloud free, this is trivial... |
---|
973 | total_albedo(:,:,:,jlev) = 0.0_jprb |
---|
974 | total_albedo(:,1,1,jlev) = total_albedo_below(:,1,1) |
---|
975 | total_albedo_direct(:,:,:,jlev) = 0.0_jprb |
---|
976 | total_albedo_direct(:,1,1,jlev) = total_albedo_below_direct(:,1,1) |
---|
977 | |
---|
978 | else if (config%i_3d_sw_entrapment == IEntrapmentMaximum & |
---|
979 | & .or. is_clear_sky_layer(jlev-1)) then |
---|
980 | ! "Maximum entrapment": use the overlap matrices u_matrix and v_matrix |
---|
981 | ! (this is the original SPARTACUS method) |
---|
982 | total_albedo(:,:,:,jlev) = singlemat_x_mat(ng,ng,nreg,& |
---|
983 | & u_matrix(:,:,jlev,jcol), & |
---|
984 | & mat_x_singlemat(ng,ng,nreg,total_albedo_below,& |
---|
985 | & v_matrix(:,:,jlev,jcol))) |
---|
986 | total_albedo_direct(:,:,:,jlev) = singlemat_x_mat(ng,ng,nreg,& |
---|
987 | & u_matrix(:,:,jlev,jcol), & |
---|
988 | & mat_x_singlemat(ng,ng,nreg,total_albedo_below_direct,& |
---|
989 | & v_matrix(:,:,jlev,jcol))) |
---|
990 | |
---|
991 | else if (config%i_3d_sw_entrapment == IEntrapmentZero) then |
---|
992 | ! "Zero entrapment": even radiation transported |
---|
993 | ! laterally between regions in the layers below is |
---|
994 | ! reflected back up into the same region. First diffuse |
---|
995 | ! radiation: |
---|
996 | total_albedo(:,:,:,jlev) = 0.0_jprb |
---|
997 | do jreg = 1,nreg ! Target layer (jlev-1) |
---|
998 | do jreg2 = 1,nreg ! Current layer (jlev) |
---|
999 | total_albedo(:,jreg,jreg,jlev) = total_albedo(:,jreg,jreg,jlev) & |
---|
1000 | & + sum(total_albedo_below(:,:,jreg2),2) & |
---|
1001 | & * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1002 | end do |
---|
1003 | end do |
---|
1004 | ! ...then direct radiation: |
---|
1005 | total_albedo_direct(:,:,:,jlev) = 0.0_jprb |
---|
1006 | do jreg = 1,nreg ! Target layer (jlev-1) |
---|
1007 | do jreg2 = 1,nreg ! Current layer (jlev) |
---|
1008 | total_albedo_direct(:,jreg,jreg,jlev) = total_albedo_direct(:,jreg,jreg,jlev) & |
---|
1009 | & + sum(total_albedo_below_direct(:,:,jreg2),2) & |
---|
1010 | & * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1011 | end do |
---|
1012 | end do |
---|
1013 | |
---|
1014 | else |
---|
1015 | ! Controlled entrapment |
---|
1016 | |
---|
1017 | #ifdef EXPLICIT_EDGE_ENTRAPMENT |
---|
1018 | ! If "EXPLICIT_EDGE_ENTRAPMENT" is defined then we use the |
---|
1019 | ! explicit entrapment approach for both horizontal transport |
---|
1020 | ! within regions, and horizontal transport between regions |
---|
1021 | ! (otherwise, horizontal transport between regions is |
---|
1022 | ! automatically treated using maximum entrapment). This is |
---|
1023 | ! experimental, which is why it is not a run-time option. |
---|
1024 | |
---|
1025 | if (config%i_3d_sw_entrapment == IEntrapmentEdgeOnly) then |
---|
1026 | #endif |
---|
1027 | ! Add the contribution from off-diagonal elements of the |
---|
1028 | ! albedo matrix in the lower layer, i.e. radiation that |
---|
1029 | ! flows between regions... |
---|
1030 | |
---|
1031 | ! First diffuse radiation: |
---|
1032 | albedo_part = total_albedo_below |
---|
1033 | do jreg = 1,nreg |
---|
1034 | albedo_part(:,jreg,jreg) = 0.0_jprb |
---|
1035 | end do |
---|
1036 | total_albedo(:,:,:,jlev) = singlemat_x_mat(ng,ng,nreg,& |
---|
1037 | & u_matrix(:,:,jlev,jcol), & |
---|
1038 | & mat_x_singlemat(ng,ng,nreg,albedo_part,& |
---|
1039 | & v_matrix(:,:,jlev,jcol))) |
---|
1040 | ! ...then direct radiation: |
---|
1041 | albedo_part = total_albedo_below_direct |
---|
1042 | do jreg = 1,nreg |
---|
1043 | albedo_part(:,jreg,jreg) = 0.0_jprb |
---|
1044 | end do |
---|
1045 | total_albedo_direct(:,:,:,jlev) = singlemat_x_mat(ng,ng,nreg,& |
---|
1046 | & u_matrix(:,:,jlev,jcol), & |
---|
1047 | & mat_x_singlemat(ng,ng,nreg,albedo_part,& |
---|
1048 | & v_matrix(:,:,jlev,jcol))) |
---|
1049 | |
---|
1050 | #ifdef EXPLICIT_EDGE_ENTRAPMENT |
---|
1051 | end if |
---|
1052 | #endif |
---|
1053 | |
---|
1054 | ! Now the contribution from the diagonals of the albedo |
---|
1055 | ! matrix in the lower layer |
---|
1056 | if (config%i_3d_sw_entrapment == IEntrapmentEdgeOnly & |
---|
1057 | & .or. (.not. config%do_3d_effects)) then |
---|
1058 | ! "Edge-only entrapment": the operation we perform is |
---|
1059 | ! essentially diag(total_albedo) += matmul(transpose(v_matrix), |
---|
1060 | ! diag(total_albedo_below)). |
---|
1061 | do jreg = 1,nreg |
---|
1062 | do jreg2 = 1,nreg |
---|
1063 | total_albedo(:,jreg,jreg,jlev) & |
---|
1064 | & = total_albedo(:,jreg,jreg,jlev) & |
---|
1065 | & + total_albedo_below(:,jreg2,jreg2) & |
---|
1066 | & * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1067 | total_albedo_direct(:,jreg,jreg,jlev) & |
---|
1068 | & = total_albedo_direct(:,jreg,jreg,jlev) & |
---|
1069 | & + total_albedo_below_direct(:,jreg2,jreg2) & |
---|
1070 | & * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1071 | end do |
---|
1072 | end do |
---|
1073 | |
---|
1074 | else |
---|
1075 | ! "Explicit entrapment" |
---|
1076 | |
---|
1077 | do jreg2 = 1,nreg |
---|
1078 | ! Loop through each region in the lower layer. For one |
---|
1079 | ! of the regions in the lower layer, we are imagining it |
---|
1080 | ! to be divided into "nreg" subregions that map onto the |
---|
1081 | ! regions in the upper layer. The rate of exchange |
---|
1082 | ! between these subregions is computed via a coupled |
---|
1083 | ! differential equation written in terms of a singular |
---|
1084 | ! exchange matrix (there are only terms for the exchange |
---|
1085 | ! between subregions, but no propagation effects since |
---|
1086 | ! we already know the albedo of this region). This is |
---|
1087 | ! solved using the matrix-exponential method. |
---|
1088 | |
---|
1089 | ! Use the following array for the transfer of either |
---|
1090 | ! diffuse or direct radiation (despite the name), per |
---|
1091 | ! unit horizontal distance travelled |
---|
1092 | transfer_rate_diffuse = 0.0_jprb |
---|
1093 | |
---|
1094 | ! As we need to reference the layer above the interface, |
---|
1095 | ! don't do the following on the highest layer |
---|
1096 | if (jlev > 1) then |
---|
1097 | |
---|
1098 | ! Given a horizontal migration distance, there is |
---|
1099 | ! still uncertainty about how much entrapment occurs |
---|
1100 | ! associated with how one assumes cloud boundaries |
---|
1101 | ! line up in adjacent layers. "overhang_factor" |
---|
1102 | ! can be varied between 0.0 (the boundaries line up to |
---|
1103 | ! the greatest extent possible given the overlap |
---|
1104 | ! parameter) and 1.0 (the boundaries line up to the |
---|
1105 | ! minimum extent possible); here this is used to |
---|
1106 | ! produce a scaling factor for the transfer rate. |
---|
1107 | transfer_scaling = 1.0_jprb - (1.0_jprb - config%overhang_factor) & |
---|
1108 | & * cloud%overlap_param(jcol,jlev-1) & |
---|
1109 | & * min(region_fracs(jreg2,jlev,jcol),region_fracs(jreg2,jlev-1,jcol)) & |
---|
1110 | & / max(config%cloud_fraction_threshold, region_fracs(jreg2,jlev,jcol)) |
---|
1111 | |
---|
1112 | do jreg = 1, nreg-1 |
---|
1113 | ! Compute lateral transfer rates from region jreg to |
---|
1114 | ! jreg+1 as before, but without length scale which |
---|
1115 | ! is wavelength dependent. |
---|
1116 | |
---|
1117 | ! Recall that overlap indexing is |
---|
1118 | ! u_matrix(upper_region, lower_region, level, |
---|
1119 | ! column). |
---|
1120 | transfer_rate_diffuse(jreg,jreg+1) = transfer_scaling & |
---|
1121 | & * edge_length(jreg,jlev-1) / max(u_matrix(jreg,jreg2,jlev,jcol),1.0e-5_jprb) |
---|
1122 | ! Compute transfer rates from region jreg+1 to jreg |
---|
1123 | transfer_rate_diffuse(jreg+1,jreg) = transfer_scaling & |
---|
1124 | & * edge_length(jreg,jlev-1) / max(u_matrix(jreg+1,jreg2,jlev,jcol),1.0e-5_jprb) |
---|
1125 | end do |
---|
1126 | |
---|
1127 | ! Compute transfer rates directly between regions 1 |
---|
1128 | ! and 3 (not used below) |
---|
1129 | if (edge_length(3,jlev) > 0.0_jprb) then |
---|
1130 | transfer_rate_diffuse(1,3) = transfer_scaling & |
---|
1131 | & * edge_length(3,jlev-1) / max(u_matrix(1,jreg2,jlev,jcol),1.0e-5_jprb) |
---|
1132 | transfer_rate_diffuse(3,1) = transfer_scaling & |
---|
1133 | & * edge_length(3,jlev-1) / max(u_matrix(3,jreg2,jlev,jcol),1.0e-5_jprb) |
---|
1134 | end if |
---|
1135 | end if |
---|
1136 | |
---|
1137 | ! Compute matrix of exchange coefficients |
---|
1138 | entrapment = 0.0_jprb |
---|
1139 | inv_effective_size = min(cloud%inv_cloud_effective_size(jcol,jlev-1), & |
---|
1140 | & 1.0_jprb/config%min_cloud_effective_size) |
---|
1141 | do jreg = 1,nreg-1 |
---|
1142 | ! Diffuse transport down and up with one random |
---|
1143 | ! scattering event |
---|
1144 | if (config%i_3d_sw_entrapment == IEntrapmentExplicit) then |
---|
1145 | fractal_factor = 1.0_jprb / sqrt(max(1.0_jprb, 2.5_jprb*x_diffuse(:,jreg2) & |
---|
1146 | & * inv_effective_size)) |
---|
1147 | entrapment(:,jreg+1,jreg) = entrapment(:,jreg+1,jreg) & |
---|
1148 | & + transfer_rate_diffuse(jreg,jreg+1)*x_diffuse(:,jreg2) & |
---|
1149 | & * fractal_factor |
---|
1150 | entrapment(:,jreg,jreg+1) = entrapment(:,jreg,jreg+1) & |
---|
1151 | & + transfer_rate_diffuse(jreg+1,jreg)*x_diffuse(:,jreg2) & |
---|
1152 | & * fractal_factor |
---|
1153 | else |
---|
1154 | entrapment(:,jreg+1,jreg) = entrapment(:,jreg+1,jreg) & |
---|
1155 | & + transfer_rate_diffuse(jreg,jreg+1)*x_diffuse(:,jreg2) |
---|
1156 | entrapment(:,jreg,jreg+1) = entrapment(:,jreg,jreg+1) & |
---|
1157 | & + transfer_rate_diffuse(jreg+1,jreg)*x_diffuse(:,jreg2) |
---|
1158 | end if |
---|
1159 | entrapment(:,jreg,jreg) = entrapment(:,jreg,jreg) & |
---|
1160 | & - entrapment(:,jreg+1,jreg) |
---|
1161 | entrapment(:,jreg+1,jreg+1) = entrapment(:,jreg+1,jreg+1) & |
---|
1162 | & - entrapment(:,jreg,jreg+1) |
---|
1163 | end do |
---|
1164 | |
---|
1165 | ! If rate of exchange is excessive the expm can throw a |
---|
1166 | ! floating point exception, even if it tends towards a |
---|
1167 | ! trival limit, so we cap the maximum input to expm by |
---|
1168 | ! scaling down if necessary |
---|
1169 | do jg = 1,ng |
---|
1170 | max_entr = -min(entrapment(jg,1,1),entrapment(jg,2,2)) |
---|
1171 | if (max_entr > config%max_cloud_od) then |
---|
1172 | ! Scale down all inputs for this g point |
---|
1173 | entrapment(jg,:,:) = entrapment(jg,:,:) * (config%max_cloud_od/max_entr) |
---|
1174 | end if |
---|
1175 | end do |
---|
1176 | |
---|
1177 | ! Since the matrix to be exponentiated has a simple |
---|
1178 | ! structure we may use a faster method described in the |
---|
1179 | ! appendix of Hogan et al. (GMD 2018) |
---|
1180 | #define USE_FAST_EXPM_EXCHANGE 1 |
---|
1181 | #ifdef USE_FAST_EXPM_EXCHANGE |
---|
1182 | if (nreg == 2) then |
---|
1183 | call fast_expm_exchange(ng, ng, entrapment(:,2,1), entrapment(:,1,2), & |
---|
1184 | & albedo_part) |
---|
1185 | else |
---|
1186 | call fast_expm_exchange(ng, ng, entrapment(:,2,1), entrapment(:,1,2), & |
---|
1187 | & entrapment(:,3,2), entrapment(:,2,3), & |
---|
1188 | & albedo_part) |
---|
1189 | end if |
---|
1190 | #else |
---|
1191 | ! Use matrix exponential to compute rate of exchange |
---|
1192 | albedo_part = entrapment |
---|
1193 | call expm(ng, ng, nreg, albedo_part, IMatrixPatternDense) |
---|
1194 | n_calls_expm = n_calls_expm + ng |
---|
1195 | #endif |
---|
1196 | |
---|
1197 | #ifndef EXPLICIT_EDGE_ENTRAPMENT |
---|
1198 | ! Scale to get the contribution to the diffuse albedo |
---|
1199 | do jreg3 = 1,nreg |
---|
1200 | do jreg = 1,nreg |
---|
1201 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) & |
---|
1202 | & * v_matrix(jreg2,jreg,jlev,jcol) * total_albedo_below(:,jreg2,jreg2) |
---|
1203 | end do |
---|
1204 | end do |
---|
1205 | #else |
---|
1206 | ! The following is an experimental treatment that tries |
---|
1207 | ! to explicitly account for the horizontal distance |
---|
1208 | ! traveled by radiation that passes through cloud sides |
---|
1209 | entrapment = albedo_part |
---|
1210 | albedo_part = 0.0_jprb |
---|
1211 | ! Scale to get the contribution to the diffuse albedo |
---|
1212 | do jreg3 = 1,nreg ! TO upper region |
---|
1213 | do jreg = 1,nreg ! FROM upper region |
---|
1214 | transfer_scaling = 1.0_jprb - (1.0_jprb - config%overhang_factor) & |
---|
1215 | & * cloud%overlap_param(jcol,jlev-1) & |
---|
1216 | & * min(region_fracs(jreg,jlev,jcol), region_fracs(jreg,jlev,jcol)) & |
---|
1217 | & / max(config%cloud_fraction_threshold, region_fracs(jreg,jlev,jcol)) |
---|
1218 | do jreg4 = 1,nreg ! VIA first lower region (jreg2 is second lower region) |
---|
1219 | if (.not. (jreg4 == jreg .and. jreg4 /= jreg2)) then |
---|
1220 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) + entrapment(:,jreg3,jreg) & |
---|
1221 | & * v_matrix(jreg4,jreg,jlev,jcol) * total_albedo_below(:,jreg2,jreg4) |
---|
1222 | else |
---|
1223 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) & |
---|
1224 | & + v_matrix(jreg4,jreg,jlev,jcol) * total_albedo_below(:,jreg2,jreg4) & |
---|
1225 | & * (transfer_scaling * entrapment(:,jreg3,jreg) & |
---|
1226 | & +((1.0_jprb-transfer_scaling) * entrapment(:,jreg3,jreg2))) |
---|
1227 | end if |
---|
1228 | end do |
---|
1229 | end do |
---|
1230 | end do |
---|
1231 | #endif |
---|
1232 | |
---|
1233 | ! Increment diffuse albedo |
---|
1234 | total_albedo(:,:,:,jlev) = total_albedo(:,:,:,jlev) + albedo_part |
---|
1235 | |
---|
1236 | ! Now do the same for the direct albedo |
---|
1237 | entrapment = 0.0_jprb |
---|
1238 | do jreg = 1,nreg-1 |
---|
1239 | ! Direct transport down and diffuse up with one random |
---|
1240 | ! scattering event |
---|
1241 | if (config%i_3d_sw_entrapment == IEntrapmentExplicit) then |
---|
1242 | fractal_factor = 1.0_jprb / sqrt(max(1.0_jprb, 2.5_jprb*x_direct(:,jreg2) & |
---|
1243 | & * inv_effective_size)) |
---|
1244 | entrapment(:,jreg+1,jreg) = entrapment(:,jreg+1,jreg) & |
---|
1245 | & + transfer_rate_diffuse(jreg,jreg+1)*x_direct(:,jreg2) & |
---|
1246 | & * fractal_factor |
---|
1247 | entrapment(:,jreg,jreg+1) = entrapment(:,jreg,jreg+1) & |
---|
1248 | & + transfer_rate_diffuse(jreg+1,jreg)*x_direct(:,jreg2) & |
---|
1249 | & * fractal_factor |
---|
1250 | else |
---|
1251 | entrapment(:,jreg+1,jreg) = entrapment(:,jreg+1,jreg) & |
---|
1252 | & + transfer_rate_diffuse(jreg,jreg+1)*x_direct(:,jreg2) |
---|
1253 | entrapment(:,jreg,jreg+1) = entrapment(:,jreg,jreg+1) & |
---|
1254 | & + transfer_rate_diffuse(jreg+1,jreg)*x_direct(:,jreg2) |
---|
1255 | end if |
---|
1256 | entrapment(:,jreg,jreg) = entrapment(:,jreg,jreg) & |
---|
1257 | & - entrapment(:,jreg+1,jreg) |
---|
1258 | entrapment(:,jreg+1,jreg+1) = entrapment(:,jreg+1,jreg+1) & |
---|
1259 | & - entrapment(:,jreg,jreg+1) |
---|
1260 | end do |
---|
1261 | |
---|
1262 | ! If rate of exchange is excessive the expm can throw a |
---|
1263 | ! floating point exception, even if it tends towards a |
---|
1264 | ! trival limit, so we cap the maximum input to expm by |
---|
1265 | ! scaling down if necessary |
---|
1266 | do jg = 1,ng |
---|
1267 | max_entr = -min(entrapment(jg,1,1),entrapment(jg,2,2)) |
---|
1268 | if (max_entr > config%max_cloud_od) then |
---|
1269 | ! Scale down all inputs for this g point |
---|
1270 | entrapment(jg,:,:) = entrapment(jg,:,:) * (config%max_cloud_od/max_entr) |
---|
1271 | end if |
---|
1272 | end do |
---|
1273 | |
---|
1274 | |
---|
1275 | #ifdef USE_FAST_EXPM_EXCHANGE |
---|
1276 | if (nreg == 2) then |
---|
1277 | call fast_expm_exchange(ng, ng, entrapment(:,2,1), entrapment(:,1,2), & |
---|
1278 | & albedo_part) |
---|
1279 | else |
---|
1280 | call fast_expm_exchange(ng, ng, entrapment(:,2,1), entrapment(:,1,2), & |
---|
1281 | & entrapment(:,3,2), entrapment(:,2,3), & |
---|
1282 | & albedo_part) |
---|
1283 | end if |
---|
1284 | #else |
---|
1285 | albedo_part = entrapment |
---|
1286 | call expm(ng, ng, nreg, albedo_part, IMatrixPatternDense) |
---|
1287 | n_calls_expm = n_calls_expm + ng |
---|
1288 | #endif |
---|
1289 | |
---|
1290 | #ifndef EXPLICIT_EDGE_ENTRAPMENT |
---|
1291 | do jreg3 = 1,nreg |
---|
1292 | do jreg = 1,nreg |
---|
1293 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) & |
---|
1294 | & * v_matrix(jreg2,jreg,jlev,jcol) * total_albedo_below_direct(:,jreg2,jreg2) |
---|
1295 | end do |
---|
1296 | end do |
---|
1297 | #else |
---|
1298 | entrapment = albedo_part |
---|
1299 | albedo_part = 0.0_jprb |
---|
1300 | do jreg3 = 1,nreg |
---|
1301 | do jreg = 1,nreg |
---|
1302 | transfer_scaling = 1.0_jprb - (1.0_jprb - config%overhang_factor) & |
---|
1303 | & * cloud%overlap_param(jcol,jlev-1) & |
---|
1304 | & * min(region_fracs(jreg,jlev,jcol), region_fracs(jreg,jlev-1,jcol)) & |
---|
1305 | & / max(config%cloud_fraction_threshold, region_fracs(jreg,jlev,jcol)) |
---|
1306 | do jreg4 = 1,nreg |
---|
1307 | if (.not. (jreg4 == jreg .and. jreg4 /= jreg2)) then |
---|
1308 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) + entrapment(:,jreg3,jreg) & |
---|
1309 | & * v_matrix(jreg4,jreg,jlev,jcol) * total_albedo_below_direct(:,jreg2,jreg4) |
---|
1310 | else |
---|
1311 | albedo_part(:,jreg3,jreg) = albedo_part(:,jreg3,jreg) & |
---|
1312 | & + v_matrix(jreg4,jreg,jlev,jcol) * total_albedo_below_direct(:,jreg2,jreg4) & |
---|
1313 | & * (transfer_scaling * entrapment(:,jreg3,jreg) & |
---|
1314 | & +((1.0_jprb-transfer_scaling) * entrapment(:,jreg3,jreg2))) |
---|
1315 | end if |
---|
1316 | end do |
---|
1317 | end do |
---|
1318 | end do |
---|
1319 | |
---|
1320 | #endif |
---|
1321 | ! Increment direct albedo |
---|
1322 | total_albedo_direct(:,:,:,jlev) = total_albedo_direct(:,:,:,jlev) + albedo_part |
---|
1323 | |
---|
1324 | end do |
---|
1325 | |
---|
1326 | end if |
---|
1327 | end if |
---|
1328 | |
---|
1329 | if ((config%i_3d_sw_entrapment == IEntrapmentExplicitNonFractal & |
---|
1330 | & .or. config%i_3d_sw_entrapment == IEntrapmentExplicit) & |
---|
1331 | & .and. .not. (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev-1))) then |
---|
1332 | ! Horizontal migration distances are averaged when |
---|
1333 | ! applying overlap rules, so equation is |
---|
1334 | ! x_above=matmul(transpose(v_matrix),x_below) |
---|
1335 | |
---|
1336 | ! We do this into temporary arrays... |
---|
1337 | x_direct_above = 0.0_jprb |
---|
1338 | x_diffuse_above = 0.0_jprb |
---|
1339 | |
---|
1340 | nregactive = nreg |
---|
1341 | if (is_clear_sky_layer(jlev)) then |
---|
1342 | nregactive = 1 |
---|
1343 | end if |
---|
1344 | |
---|
1345 | do jreg = 1,nreg ! Target layer (jlev-1) |
---|
1346 | do jreg2 = 1,nregactive ! Current layer (jlev) |
---|
1347 | x_direct_above(:,jreg) = x_direct_above(:,jreg) & |
---|
1348 | & + x_direct(:,jreg2) * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1349 | x_diffuse_above(:,jreg) = x_diffuse_above(:,jreg) & |
---|
1350 | & + x_diffuse(:,jreg2) * v_matrix(jreg2,jreg,jlev,jcol) |
---|
1351 | end do |
---|
1352 | end do |
---|
1353 | |
---|
1354 | !... then copy out of the temporary arrays |
---|
1355 | x_direct = x_direct_above |
---|
1356 | x_diffuse = x_diffuse_above |
---|
1357 | end if |
---|
1358 | |
---|
1359 | end do ! Reverse loop over levels |
---|
1360 | |
---|
1361 | ! -------------------------------------------------------- |
---|
1362 | ! Section 5: Compute fluxes |
---|
1363 | ! -------------------------------------------------------- |
---|
1364 | |
---|
1365 | ! Top-of-atmosphere fluxes into the regions of the top-most |
---|
1366 | ! layer, zero since we assume no diffuse downwelling |
---|
1367 | flux_dn_below = 0.0_jprb |
---|
1368 | ! Direct downwelling flux (into a plane perpendicular to the |
---|
1369 | ! sun) entering the top of each region in the top-most layer |
---|
1370 | do jreg = 1,nreg |
---|
1371 | direct_dn_below(:,jreg) = incoming_sw(:,jcol)*region_fracs(jreg,1,jcol) |
---|
1372 | end do |
---|
1373 | ! We're using flux_up_above as a container; actually its |
---|
1374 | ! interpretation at top of atmosphere here is just 'below' the |
---|
1375 | ! TOA interface, so using the regions of the first model layer |
---|
1376 | flux_up_above = mat_x_vec(ng,ng,nreg,total_albedo_direct(:,:,:,1),direct_dn_below) |
---|
1377 | |
---|
1378 | if (config%do_clear) then |
---|
1379 | flux_dn_clear = 0.0_jprb |
---|
1380 | direct_dn_clear(:) = incoming_sw(:,jcol) |
---|
1381 | flux_up_clear = direct_dn_clear*total_albedo_clear_direct(:,1) |
---|
1382 | end if |
---|
1383 | |
---|
1384 | ! Store the TOA broadband fluxes |
---|
1385 | flux%sw_up(jcol,1) = sum(sum(flux_up_above,1)) |
---|
1386 | flux%sw_dn(jcol,1) = mu0 * sum(direct_dn_clear(:)) |
---|
1387 | if (allocated(flux%sw_dn_direct)) then |
---|
1388 | flux%sw_dn_direct(jcol,1) = flux%sw_dn(jcol,1) |
---|
1389 | end if |
---|
1390 | if (config%do_clear) then |
---|
1391 | flux%sw_up_clear(jcol,1) = sum(flux_up_clear) |
---|
1392 | flux%sw_dn_clear(jcol,1) = flux%sw_dn(jcol,1) |
---|
1393 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
1394 | flux%sw_dn_direct_clear(jcol,1) = flux%sw_dn_clear(jcol,1) |
---|
1395 | end if |
---|
1396 | end if |
---|
1397 | |
---|
1398 | ! Save the spectral fluxes if required |
---|
1399 | if (config%do_save_spectral_flux) then |
---|
1400 | call indexed_sum(sum(flux_up_above(:,:),2), & |
---|
1401 | & config%i_spec_from_reordered_g_sw, & |
---|
1402 | & flux%sw_up_band(:,jcol,1)) |
---|
1403 | call indexed_sum(sum(direct_dn_below(:,:),2), & |
---|
1404 | & config%i_spec_from_reordered_g_sw, & |
---|
1405 | & flux%sw_dn_band(:,jcol,1)) |
---|
1406 | flux%sw_dn_band(:,jcol,1) = mu0 * flux%sw_dn_band(:,jcol,1) |
---|
1407 | if (allocated(flux%sw_dn_direct_band)) then |
---|
1408 | flux%sw_dn_direct_band(:,jcol,1) = flux%sw_dn_band(:,jcol,1) |
---|
1409 | end if |
---|
1410 | if (config%do_clear) then |
---|
1411 | flux%sw_dn_clear_band(:,jcol,1) = flux%sw_dn_band(:,jcol,1) |
---|
1412 | call indexed_sum(flux_up_clear, & |
---|
1413 | & config%i_spec_from_reordered_g_sw, & |
---|
1414 | & flux%sw_up_clear_band(:,jcol,1)) |
---|
1415 | if (allocated(flux%sw_dn_direct_clear_band)) then |
---|
1416 | flux%sw_dn_direct_clear_band(:,jcol,1) & |
---|
1417 | & = flux%sw_dn_clear_band(:,jcol,1) |
---|
1418 | end if |
---|
1419 | end if |
---|
1420 | end if |
---|
1421 | |
---|
1422 | ! Final loop back down through the atmosphere to compute fluxes |
---|
1423 | do jlev = 1,nlev |
---|
1424 | |
---|
1425 | #ifdef PRINT_ENTRAPMENT_DATA |
---|
1426 | if (config%i_3d_sw_entrapment == IEntrapmentExplicitNonFractal & |
---|
1427 | & .or. config%i_3d_sw_entrapment == IEntrapmentExplicit) then |
---|
1428 | ! Save downwelling direct and diffuse fluxes at the top of |
---|
1429 | ! layer "jlev" in each of the regions of layer "jlev" |
---|
1430 | if (nreg == 2) then |
---|
1431 | write(102,'(i4,i4,4e14.6)') jcol, jlev, direct_dn_below(1,:), flux_dn_below(1,:) |
---|
1432 | else |
---|
1433 | write(102,'(i4,i4,6e14.6)') jcol, jlev, direct_dn_below(1,1:3), flux_dn_below(1,1:3) |
---|
1434 | end if |
---|
1435 | end if |
---|
1436 | #endif |
---|
1437 | |
---|
1438 | ! Compute the solar downwelling "source" at the base of the |
---|
1439 | ! layer due to scattering of the direct beam within it |
---|
1440 | if (config%do_clear) then |
---|
1441 | source_dn_clear = trans_dir_diff_clear(:,jlev)*direct_dn_clear |
---|
1442 | end if |
---|
1443 | source_dn(:,:) = mat_x_vec(ng,ng,nreg,trans_dir_diff(:,:,:,jlev),direct_dn_below, & |
---|
1444 | & is_clear_sky_layer(jlev)) |
---|
1445 | |
---|
1446 | ! Compute direct downwelling flux in each region at base of |
---|
1447 | ! current layer |
---|
1448 | if (config%do_clear) then |
---|
1449 | direct_dn_clear = trans_dir_dir_clear(:,jlev)*direct_dn_clear |
---|
1450 | end if |
---|
1451 | direct_dn_above = mat_x_vec(ng,ng,nreg,trans_dir_dir(:,:,:,jlev),direct_dn_below, & |
---|
1452 | & is_clear_sky_layer(jlev)) |
---|
1453 | |
---|
1454 | ! Integrate downwelling direct flux across spectrum and |
---|
1455 | ! regions, and store (the diffuse part will be added later) |
---|
1456 | flux%sw_dn(jcol,jlev+1) = mu0 * sum(sum(direct_dn_above,1)) |
---|
1457 | if (allocated(flux%sw_dn_direct)) then |
---|
1458 | flux%sw_dn_direct(jcol,jlev+1) = flux%sw_dn(jcol,jlev+1) |
---|
1459 | end if |
---|
1460 | if (config%do_clear) then |
---|
1461 | flux%sw_dn_clear(jcol,jlev+1) = mu0 * sum(direct_dn_clear) |
---|
1462 | if (allocated(flux%sw_dn_direct_clear)) then |
---|
1463 | flux%sw_dn_direct_clear(jcol,jlev+1) & |
---|
1464 | & = flux%sw_dn_clear(jcol,jlev+1) |
---|
1465 | end if |
---|
1466 | end if |
---|
1467 | |
---|
1468 | if (config%do_save_spectral_flux) then |
---|
1469 | call indexed_sum(sum(direct_dn_above,2), & |
---|
1470 | & config%i_spec_from_reordered_g_sw, & |
---|
1471 | & flux%sw_dn_band(:,jcol,jlev+1)) |
---|
1472 | flux%sw_dn_band(:,jcol,jlev+1) = mu0 * flux%sw_dn_band(:,jcol,jlev+1) |
---|
1473 | |
---|
1474 | if (allocated(flux%sw_dn_direct_band)) then |
---|
1475 | flux%sw_dn_direct_band(:,jcol,jlev+1) & |
---|
1476 | & = flux%sw_dn_band(:,jcol,jlev+1) |
---|
1477 | end if |
---|
1478 | if (config%do_clear) then |
---|
1479 | call indexed_sum(direct_dn_clear, & |
---|
1480 | & config%i_spec_from_reordered_g_sw, & |
---|
1481 | & flux%sw_dn_clear_band(:,jcol,jlev+1)) |
---|
1482 | flux%sw_dn_clear_band(:,jcol,jlev+1) = mu0 & |
---|
1483 | & * flux%sw_dn_clear_band(:,jcol,jlev+1) |
---|
1484 | if (allocated(flux%sw_dn_direct_clear_band)) then |
---|
1485 | flux%sw_dn_direct_clear_band(:,jcol,jlev+1) & |
---|
1486 | & = flux%sw_dn_clear_band(:,jcol,jlev+1) |
---|
1487 | end if |
---|
1488 | end if |
---|
1489 | end if |
---|
1490 | |
---|
1491 | if (config%do_clear) then |
---|
1492 | ! Scalar operations for clear-sky fluxes |
---|
1493 | flux_dn_clear(:) = (trans_clear(:,jlev)*flux_dn_clear(:) & |
---|
1494 | & + ref_clear(:,jlev)*total_albedo_clear_direct(:,jlev+1)*direct_dn_clear & |
---|
1495 | & + source_dn_clear) & |
---|
1496 | & / (1.0_jprb - ref_clear(:,jlev)*total_albedo_clear(:,jlev+1)) |
---|
1497 | flux_up_clear(:) = total_albedo_clear_direct(:,jlev+1)*direct_dn_clear & |
---|
1498 | & + total_albedo_clear(:,jlev+1)*flux_dn_clear |
---|
1499 | end if |
---|
1500 | |
---|
1501 | if (is_clear_sky_layer(jlev)) then |
---|
1502 | ! Scalar operations for clear-sky layer |
---|
1503 | flux_dn_above(:,1) = (transmittance(:,1,1,jlev)*flux_dn_below(:,1) & |
---|
1504 | & + reflectance(:,1,1,jlev)*total_albedo_direct(:,1,1,jlev+1)*direct_dn_above(:,1) & |
---|
1505 | & + source_dn(:,1)) & |
---|
1506 | & / (1.0_jprb - reflectance(:,1,1,jlev)*total_albedo(:,1,1,jlev+1)) |
---|
1507 | flux_dn_above(:,2:nreg) = 0.0_jprb |
---|
1508 | flux_up_above(:,1) = total_albedo_direct(:,1,1,jlev+1)*direct_dn_above(:,1) & |
---|
1509 | & + total_albedo(:,1,1,jlev+1)*flux_dn_above(:,1) |
---|
1510 | flux_up_above(:,2:nreg) = 0.0_jprb |
---|
1511 | else |
---|
1512 | ! Matrix operations for cloudy layer |
---|
1513 | denominator = identity_minus_mat_x_mat(ng,ng,nreg,reflectance(:,:,:,jlev), & |
---|
1514 | & total_albedo(:,:,:,jlev+1)) |
---|
1515 | total_source = mat_x_vec(ng,ng,nreg,total_albedo_direct(:,:,:,jlev+1),direct_dn_above) |
---|
1516 | |
---|
1517 | flux_dn_above = solve_vec(ng,ng,nreg,denominator, & |
---|
1518 | & mat_x_vec(ng,ng,nreg,transmittance(:,:,:,jlev),flux_dn_below) & |
---|
1519 | & + mat_x_vec(ng,ng,nreg,reflectance(:,:,:,jlev), total_source(:,:)) & |
---|
1520 | & + source_dn(:,:)) |
---|
1521 | flux_up_above = mat_x_vec(ng,ng,nreg,total_albedo(:,:,:,jlev+1), & |
---|
1522 | & flux_dn_above) + total_source(:,:) |
---|
1523 | end if |
---|
1524 | |
---|
1525 | ! Account for overlap rules in translating fluxes just above |
---|
1526 | ! a layer interface to the values just below |
---|
1527 | if (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev+1)) then |
---|
1528 | ! Regions in current layer map directly on to regions in |
---|
1529 | ! layer below |
---|
1530 | flux_dn_below = flux_dn_above |
---|
1531 | direct_dn_below = direct_dn_above |
---|
1532 | else |
---|
1533 | ! Apply downward overlap matrix to compute direct |
---|
1534 | ! downwelling flux entering the top of each region in the |
---|
1535 | ! layer below |
---|
1536 | flux_dn_below = singlemat_x_vec(ng,ng,nreg,v_matrix(:,:,jlev+1,jcol), & |
---|
1537 | & flux_dn_above) |
---|
1538 | direct_dn_below = singlemat_x_vec(ng,ng,nreg,v_matrix(:,:,jlev+1,jcol), & |
---|
1539 | & direct_dn_above) |
---|
1540 | end if |
---|
1541 | |
---|
1542 | ! Store the broadband fluxes |
---|
1543 | flux%sw_up(jcol,jlev+1) = sum(sum(flux_up_above,1)) |
---|
1544 | flux%sw_dn(jcol,jlev+1) & |
---|
1545 | & = flux%sw_dn(jcol,jlev+1) + sum(sum(flux_dn_above,1)) |
---|
1546 | if (config%do_clear) then |
---|
1547 | flux%sw_up_clear(jcol,jlev+1) = sum(flux_up_clear) |
---|
1548 | flux%sw_dn_clear(jcol,jlev+1) & |
---|
1549 | & = flux%sw_dn_clear(jcol,jlev+1) + sum(flux_dn_clear) |
---|
1550 | end if |
---|
1551 | |
---|
1552 | ! Save the spectral fluxes if required |
---|
1553 | if (config%do_save_spectral_flux) then |
---|
1554 | call indexed_sum(sum(flux_up_above,2), & |
---|
1555 | & config%i_spec_from_reordered_g_sw, & |
---|
1556 | & flux%sw_up_band(:,jcol,jlev+1)) |
---|
1557 | call add_indexed_sum(sum(flux_dn_above,2), & |
---|
1558 | & config%i_spec_from_reordered_g_sw, & |
---|
1559 | & flux%sw_dn_band(:,jcol,jlev+1)) |
---|
1560 | if (config%do_clear) then |
---|
1561 | call indexed_sum(flux_up_clear, & |
---|
1562 | & config%i_spec_from_reordered_g_sw, & |
---|
1563 | & flux%sw_up_clear_band(:,jcol,jlev+1)) |
---|
1564 | call add_indexed_sum(flux_dn_clear, & |
---|
1565 | & config%i_spec_from_reordered_g_sw, & |
---|
1566 | & flux%sw_dn_clear_band(:,jcol,jlev+1)) |
---|
1567 | end if |
---|
1568 | end if |
---|
1569 | |
---|
1570 | end do ! Final loop over levels |
---|
1571 | |
---|
1572 | ! Store surface spectral fluxes, if required (after the end of |
---|
1573 | ! the final loop over levels, the current values of these arrays |
---|
1574 | ! will be the surface values) |
---|
1575 | flux%sw_dn_diffuse_surf_g(:,jcol) = sum(flux_dn_above,2) |
---|
1576 | flux%sw_dn_direct_surf_g(:,jcol) = mu0 * sum(direct_dn_above,2) |
---|
1577 | if (config%do_clear) then |
---|
1578 | flux%sw_dn_diffuse_surf_clear_g(:,jcol) = flux_dn_clear |
---|
1579 | flux%sw_dn_direct_surf_clear_g(:,jcol) = mu0 * direct_dn_clear |
---|
1580 | end if |
---|
1581 | |
---|
1582 | end do ! Loop over columns |
---|
1583 | if (config%iverbose >= 3) then |
---|
1584 | write(nulout,*) |
---|
1585 | end if |
---|
1586 | |
---|
1587 | ! Report number of calls to each method of solving single-layer |
---|
1588 | ! two-stream equations |
---|
1589 | if (config%iverbose >= 4) then |
---|
1590 | write(nulout,'(a,i0)') ' Matrix-exponential calls: ', n_calls_expm |
---|
1591 | write(nulout,'(a,i0)') ' Meador-Weaver calls: ', n_calls_meador_weaver |
---|
1592 | end if |
---|
1593 | |
---|
1594 | if (lhook) call dr_hook('radiation_spartacus_sw:solver_spartacus_sw',1,hook_handle) |
---|
1595 | |
---|
1596 | end subroutine solver_spartacus_sw |
---|
1597 | |
---|
1598 | |
---|
1599 | ! Step the horizontal migration distances from the base of a layer |
---|
1600 | ! to the top, accounting for the extra distance travelled within the |
---|
1601 | ! layer |
---|
1602 | subroutine step_migrations(ng, nreg, cloud_frac, & |
---|
1603 | & layer_depth, tan_diffuse_angle_3d, tan_sza, & |
---|
1604 | & reflectance, transmittance, ref_dir, trans_dir_dir, & |
---|
1605 | & trans_dir_diff, total_albedo_diff, total_albedo_dir, & |
---|
1606 | & x_diffuse, x_direct) |
---|
1607 | |
---|
1608 | use parkind1, only : jprb |
---|
1609 | |
---|
1610 | implicit none |
---|
1611 | |
---|
1612 | ! Inputs |
---|
1613 | |
---|
1614 | ! Number of g points and regions |
---|
1615 | integer, intent(in) :: ng, nreg |
---|
1616 | ! Cloud fraction |
---|
1617 | real(jprb), intent(in) :: cloud_frac |
---|
1618 | ! Layer depth (m), tangent of diffuse zenith angle and tangent of |
---|
1619 | ! solar zenith angle |
---|
1620 | real(jprb), intent(in) :: layer_depth, tan_diffuse_angle_3d, tan_sza |
---|
1621 | ! Reflectance and transmittance to diffuse downwelling radiation |
---|
1622 | real(jprb), intent(in), dimension(ng, nreg, nreg) :: reflectance, transmittance |
---|
1623 | ! Reflectance and transmittance to direct downwelling radiation |
---|
1624 | real(jprb), intent(in), dimension(ng, nreg, nreg) :: ref_dir, trans_dir_dir |
---|
1625 | ! Transmittance involving direct entering a layer from the top and |
---|
1626 | ! diffuse leaving from the bottom |
---|
1627 | real(jprb), intent(in), dimension(ng, nreg, nreg) :: trans_dir_diff |
---|
1628 | |
---|
1629 | ! Total albedo of direct and diffuse radiation of the atmosphere |
---|
1630 | ! below the layer in question |
---|
1631 | real(jprb), intent(in), dimension(ng, nreg, nreg) & |
---|
1632 | & :: total_albedo_diff, total_albedo_dir |
---|
1633 | |
---|
1634 | ! Inputs/outputs |
---|
1635 | |
---|
1636 | ! Horizontal migration distance (m) of reflected light |
---|
1637 | real(jprb), intent(inout), dimension(ng, nreg) :: x_diffuse, x_direct |
---|
1638 | |
---|
1639 | ! Local variables |
---|
1640 | |
---|
1641 | ! Top albedo, i.e. the albedo of the top of the layer assuming no |
---|
1642 | ! lateral transport |
---|
1643 | real(jprb), dimension(ng) :: top_albedo |
---|
1644 | |
---|
1645 | ! Multiple-scattering amplitude enhancement |
---|
1646 | real(jprb), dimension(ng) :: ms_enhancement |
---|
1647 | |
---|
1648 | ! Multiple-scattering distance enhancement |
---|
1649 | real(jprb), dimension(ng) :: x_enhancement |
---|
1650 | |
---|
1651 | |
---|
1652 | real(jprb) :: x_layer_diffuse, x_layer_direct |
---|
1653 | integer :: jreg, istartreg, iendreg |
---|
1654 | |
---|
1655 | istartreg = 1 |
---|
1656 | iendreg = nreg |
---|
1657 | |
---|
1658 | if (cloud_frac <= 0.0_jprb) then |
---|
1659 | ! Clear-sky layer: don't waste time on cloudy regions |
---|
1660 | iendreg = 1 |
---|
1661 | else if (cloud_frac >= 1.0_jprb) then |
---|
1662 | ! Overcast layer: don't waste time on clear region |
---|
1663 | istartreg = 2 |
---|
1664 | end if |
---|
1665 | |
---|
1666 | ! This is the mean horizontal distance travelled by diffuse |
---|
1667 | ! radiation that travels from the top of a layer to the centre and |
---|
1668 | ! is then scattered back up and out |
---|
1669 | x_layer_diffuse = layer_depth * tan_diffuse_angle_3d/sqrt(2.0_jprb) |
---|
1670 | |
---|
1671 | ! This is the mean horizontal distance travelled by direct |
---|
1672 | ! radiation that travels from the top of a layer to the centre and |
---|
1673 | ! is then scattered back up and out |
---|
1674 | x_layer_direct = layer_depth * sqrt(tan_sza*tan_sza & |
---|
1675 | & + tan_diffuse_angle_3d*tan_diffuse_angle_3d) * 0.5_jprb |
---|
1676 | |
---|
1677 | do jreg = istartreg,iendreg |
---|
1678 | ! Geometric series enhancement due to multiple scattering: the |
---|
1679 | ! amplitude enhancement is equal to the limit of |
---|
1680 | ! T*[1+RA+(RA)^2+(RA)^3+...] |
---|
1681 | ms_enhancement = transmittance(:,jreg,jreg) & |
---|
1682 | & / (1.0_jprb - reflectance(:,jreg,jreg)*total_albedo_diff(:,jreg,jreg)) |
---|
1683 | ! ...and the distance enhancement is approximately equal to the |
---|
1684 | ! limit of T*[1+sqrt(2)*RA+sqrt(3)*(RA)^2+sqrt(4)*(RA)^3+...] |
---|
1685 | x_enhancement = (1.0_jprb - reflectance(:,jreg,jreg)*total_albedo_diff(:,jreg,jreg))**(-1.5_jprb) |
---|
1686 | |
---|
1687 | ! Horizontal migration of direct downwelling radiation |
---|
1688 | top_albedo = max(1.0e-8_jprb, ref_dir(:,jreg,jreg) + ms_enhancement & |
---|
1689 | & * (trans_dir_diff(:,jreg,jreg)*total_albedo_diff(:,jreg,jreg) & |
---|
1690 | & +trans_dir_dir(:,jreg,jreg)*total_albedo_dir(:,jreg,jreg))) |
---|
1691 | ! The following is approximate and has been found to |
---|
1692 | ! occasionally go negative |
---|
1693 | x_direct(:,jreg) = max(0.0_jprb, x_layer_direct & |
---|
1694 | & + ((trans_dir_diff(:,jreg,jreg)*total_albedo_diff(:,jreg,jreg)*x_enhancement & |
---|
1695 | & +trans_dir_dir(:,jreg,jreg)*total_albedo_dir(:,jreg,jreg)*(x_enhancement-1.0_jprb)) & |
---|
1696 | & *(x_diffuse(:,jreg)+x_layer_diffuse) & |
---|
1697 | & +trans_dir_dir(:,jreg,jreg)*total_albedo_dir(:,jreg,jreg) & |
---|
1698 | & *(x_direct(:,jreg)+x_layer_direct)) & |
---|
1699 | & * transmittance(:,jreg,jreg) / top_albedo) |
---|
1700 | |
---|
1701 | ! Horizontal migration of diffuse downwelling radiation |
---|
1702 | top_albedo = max(1.0e-8_jprb, reflectance(:,jreg,jreg) & |
---|
1703 | & + ms_enhancement*transmittance(:,jreg,jreg)*total_albedo_diff(:,jreg,jreg)) |
---|
1704 | x_diffuse(:,jreg) = x_layer_diffuse + x_enhancement*total_albedo_diff(:,jreg,jreg) & |
---|
1705 | & *(transmittance(:,jreg,jreg)*transmittance(:,jreg,jreg)) & |
---|
1706 | & * (x_diffuse(:,jreg) + x_layer_diffuse) / top_albedo |
---|
1707 | |
---|
1708 | end do |
---|
1709 | if (iendreg < nreg) then |
---|
1710 | x_diffuse(:,iendreg+1:nreg) = 0.0_jprb |
---|
1711 | x_direct(:,iendreg+1:nreg) = 0.0_jprb |
---|
1712 | else if (istartreg == 2) then |
---|
1713 | x_diffuse(:,1) = 0.0_jprb |
---|
1714 | x_direct(:,1) = 0.0_jprb |
---|
1715 | end if |
---|
1716 | |
---|
1717 | end subroutine step_migrations |
---|
1718 | |
---|
1719 | end module radiation_spartacus_sw |
---|