source: LMDZ6/branches/LMDZ-ECRAD/libf/phylmd/ecrad/radiation_spartacus_sw.F90 @ 3880

Last change on this file since 3880 was 3880, checked in by idelkadi, 3 years ago

Online implementation of the radiative transfer code ECRAD in LMDZ.

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