source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_tripleclouds_lw.F90.or @ 5418

Last change on this file since 5418 was 4946, checked in by idelkadi, 7 months ago

The addition of explicit loops with the "omp simd reduction" directive to solve the slowness problem linked to the "sum" command (svn4938 revesion) led to non-reproducibility in MPI and mixed MPI-OMP modes in the case of Tripleclouds and Mcica solvers.
We return to the svn4848 versions of the radiation_tripleclouds_*w.F90 and radiation_mcica_*w.F90 routines.

File size: 25.6 KB
Line 
1! radiation_tripleclouds_lw.F90 - Longwave "Tripleclouds" solver
2!
3! (C) Copyright 2016- 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-28  R. Hogan  Receive emission/albedo rather than planck/emissivity
17!   2017-04-22  R. Hogan  Store surface fluxes at all g-points
18!   2017-10-23  R. Hogan  Renamed single-character variables
19!   2018-10-08  R. Hogan  Call calc_region_properties
20!   2020-09-18  R. Hogan  Replaced some array expressions with loops
21!   2020-09-19  R. Hogan  Implement the cloud-only-scattering optimization
22
23module radiation_tripleclouds_lw
24
25  public
26
27contains
28  ! Small routine for scaling cloud optical depth in the cloudy
29  ! regions
30#include "radiation_optical_depth_scaling.h"
31
32  !---------------------------------------------------------------------
33  ! This module contains just one subroutine, the longwave
34  ! "Tripleclouds" solver in which cloud inhomogeneity is treated by
35  ! dividing each model level into three regions, one clear and two
36  ! cloudy (with differing optical depth). This approach was described
37  ! by Shonk and Hogan (2008).
38  subroutine solver_tripleclouds_lw(nlev,istartcol,iendcol, &
39       &  config, cloud, &
40       &  od, ssa, g, od_cloud, ssa_cloud, g_cloud, planck_hl, &
41       &  emission, albedo, &
42       &  flux)
43
44    use parkind1, only           : jprb
45    use yomhook,  only           : lhook, dr_hook, jphook
46
47!    use radiation_io, only             : nulout
48    use radiation_config, only         : config_type, IPdfShapeGamma
49    use radiation_cloud, only          : cloud_type
50    use radiation_regions, only        : calc_region_properties
51    use radiation_overlap, only        : calc_overlap_matrices
52    use radiation_flux, only           : flux_type, indexed_sum
53    use radiation_matrix, only         : singlemat_x_vec
54    use radiation_two_stream, only     : calc_ref_trans_lw, &
55         &                               calc_no_scattering_transmittance_lw
56    use radiation_adding_ica_lw, only  : adding_ica_lw, calc_fluxes_no_scattering_lw
57    use radiation_lw_derivatives, only : calc_lw_derivatives_region
58
59    implicit none
60
61    ! Inputs
62    integer, intent(in) :: nlev               ! number of model levels
63    integer, intent(in) :: istartcol, iendcol ! range of columns to process
64    type(config_type),        intent(in) :: config
65    type(cloud_type),         intent(in) :: cloud
66
67    ! Gas and aerosol optical depth of each layer at each longwave
68    ! g-point
69    real(jprb), intent(in), dimension(config%n_g_lw,nlev,istartcol:iendcol) :: od
70
71    ! Gas and aerosol single-scattering albedo and asymmetry factor,
72    ! only if longwave scattering by aerosols is to be represented
73    real(jprb), intent(in), &
74         &  dimension(config%n_g_lw_if_scattering,nlev,istartcol:iendcol) :: ssa, g
75
76    ! Cloud and precipitation optical depth of each layer in each
77    ! longwave band
78    real(jprb), intent(in) :: od_cloud(config%n_bands_lw,nlev,istartcol:iendcol)
79
80    ! Cloud and precipitation single-scattering albedo and asymmetry
81    ! factor, only if longwave scattering by clouds is to be
82    ! represented
83    real(jprb), intent(in), dimension(config%n_bands_lw_if_scattering, &
84         &                            nlev,istartcol:iendcol) :: ssa_cloud, g_cloud
85
86    ! Planck function (emitted flux from a black body) at half levels
87    ! and at the surface at each longwave g-point
88    real(jprb), intent(in), dimension(config%n_g_lw,nlev+1,istartcol:iendcol) :: planck_hl
89
90    ! Emission (Planck*emissivity) and albedo (1-emissivity) at the
91    ! surface at each longwave g-point
92    real(jprb), intent(in), dimension(config%n_g_lw, istartcol:iendcol) :: emission, albedo
93
94    ! Optical depth, single scattering albedo and asymmetry factor in
95    ! each g-point including gas, aerosol and clouds
96    real(jprb), dimension(config%n_g_lw) :: od_total, ssa_total, g_total
97
98    ! Modified optical depth after Tripleclouds scaling to represent
99    ! cloud inhomogeneity
100    real(jprb), dimension(config%n_g_lw) :: od_cloud_new
101
102    ! Output
103    type(flux_type), intent(inout):: flux
104
105    ! Local constants
106    integer, parameter :: nregions = 3
107
108    ! In a clear-sky layer this will be 1, otherwise equal to nregions
109    integer :: nreg
110
111    ! Local variables
112   
113    ! The area fractions of each region
114    real(jprb) :: region_fracs(1:nregions,nlev,istartcol:iendcol)
115
116    ! The scaling used for the optical depth in the cloudy regions
117    real(jprb) :: od_scaling(2:nregions,nlev,istartcol:iendcol)
118
119    ! Directional overlap matrices defined at all layer interfaces
120    ! including top-of-atmosphere and the surface
121    real(jprb), dimension(nregions,nregions,nlev+1, &
122         &                istartcol:iendcol) :: u_matrix, v_matrix
123
124    ! Diffuse reflection and transmission matrices of each layer
125    real(jprb), dimension(config%n_g_lw, nregions, nlev) :: reflectance, transmittance
126
127    ! Emission by a layer into the upwelling or downwelling diffuse
128    ! streams
129    real(jprb), dimension(config%n_g_lw, nregions, nlev) &
130         &  :: source_up, source_dn
131
132    ! Clear-sky reflectance and transmittance
133    real(jprb), dimension(config%n_g_lw, nlev) &
134         &  :: ref_clear, trans_clear
135
136    ! ...clear-sky equivalent
137    real(jprb), dimension(config%n_g_lw, nlev) &
138         &  :: source_up_clear, source_dn_clear
139
140    ! Total albedo of the atmosphere/surface just above a layer
141    ! interface with respect to downwelling diffuse radiation at that
142    ! interface, where level index = 1 corresponds to the
143    ! top-of-atmosphere
144    real(jprb), dimension(config%n_g_lw, nregions, nlev+1) :: total_albedo
145
146    ! Upwelling radiation just above a layer interface due to emission
147    ! below that interface, where level index = 1 corresponds to the
148    ! top-of-atmosphere
149    real(jprb), dimension(config%n_g_lw, nregions, nlev+1) :: total_source
150
151    ! Total albedo and source of the atmosphere just below a layer interface
152    real(jprb), dimension(config%n_g_lw, nregions) &
153         &  :: total_albedo_below, total_source_below
154
155    ! Downwelling flux below and above an interface between
156    ! layers into a plane perpendicular to the direction of the sun
157    real(jprb), dimension(config%n_g_lw, nregions) &
158         &  :: flux_dn, flux_dn_below, flux_up
159
160    ! ...clear-sky equivalent (no distinction between "above/below")
161    real(jprb), dimension(config%n_g_lw, nlev+1) &
162         &  :: flux_dn_clear, flux_up_clear
163
164    ! Clear-sky equivalent, but actually its reciprocal to replace
165    ! some divisions by multiplications
166    real(jprb), dimension(config%n_g_lw, nregions) :: inv_denom
167
168    ! Identify clear-sky layers, with pseudo layers for outer space
169    ! and below the ground, both treated as single-region clear skies
170    logical :: is_clear_sky_layer(0:nlev+1)
171
172    ! Temporaries to speed up summations
173    real(jprb) :: sum_dn, sum_up
174   
175    ! Index of the highest cloudy layer
176    integer :: i_cloud_top
177
178    integer :: jcol, jlev, jg, jreg, jreg2, ng
179
180    real(jphook) :: hook_handle
181
182    if (lhook) call dr_hook('radiation_tripleclouds_lw:solver_tripleclouds_lw',0,hook_handle)
183
184    ! --------------------------------------------------------
185    ! Section 1: Prepare general variables and arrays
186    ! --------------------------------------------------------
187    ! Copy array dimensions to local variables for convenience
188    ng   = config%n_g_lw
189
190    ! Compute the wavelength-independent region fractions and
191    ! optical-depth scalings
192    call calc_region_properties(nlev,nregions,istartcol,iendcol, &
193         &  config%i_cloud_pdf_shape == IPdfShapeGamma, &
194         &  cloud%fraction, cloud%fractional_std, region_fracs, &
195         &  od_scaling, config%cloud_fraction_threshold)
196
197    ! Compute wavelength-independent overlap matrices u_matrix and v_matrix
198    call calc_overlap_matrices(nlev,nregions,istartcol,iendcol, &
199         &  region_fracs, cloud%overlap_param, &
200         &  u_matrix, v_matrix, &
201         &  decorrelation_scaling=config%cloud_inhom_decorr_scaling, &
202         &  cloud_fraction_threshold=config%cloud_fraction_threshold, &
203         &  use_beta_overlap=config%use_beta_overlap, &
204         &  cloud_cover=flux%cloud_cover_lw)
205
206    ! Main loop over columns
207    do jcol = istartcol, iendcol
208      ! --------------------------------------------------------
209      ! Section 2: Prepare column-specific variables and arrays
210      ! --------------------------------------------------------
211
212      ! Define which layers contain cloud; assume that
213      ! cloud%crop_cloud_fraction has already been called
214      is_clear_sky_layer = .true.
215      i_cloud_top = nlev+1
216      do jlev = 1,nlev
217        if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
218          is_clear_sky_layer(jlev) = .false.
219          ! Get index to the first cloudy layer from the top
220          if (i_cloud_top > jlev) then
221            i_cloud_top = jlev
222          end if
223        end if
224      end do
225      if (config%do_lw_aerosol_scattering) then
226        ! This is actually the first layer in which we need to
227        ! consider scattering
228        i_cloud_top = 1
229      end if
230
231      ! --------------------------------------------------------
232      ! Section 3: Clear-sky calculation
233      ! --------------------------------------------------------
234
235      if (.not. config%do_lw_aerosol_scattering) then
236        ! No scattering in clear-sky flux calculation; note that here
237        ! the first two dimensions of the input arrays are unpacked
238        ! into vectors inside the routine       
239        call calc_no_scattering_transmittance_lw(ng*nlev, od(:,:,jcol), &
240             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1, jcol), &
241             &  trans_clear, source_up_clear, source_dn_clear)
242        ! Ensure that clear-sky reflectance is zero since it may be
243        ! used in cloudy-sky case
244        ref_clear = 0.0_jprb
245        ! Simple down-then-up method to compute fluxes
246        call calc_fluxes_no_scattering_lw(ng, nlev, &
247             &  trans_clear, source_up_clear, source_dn_clear, &
248             &  emission(:,jcol), albedo(:,jcol), &
249             &  flux_up_clear, flux_dn_clear)
250      else
251        ! Scattering in clear-sky flux calculation
252        call calc_ref_trans_lw(ng*nlev, &
253             &  od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), &
254             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1,jcol), &
255             &  ref_clear, trans_clear, &
256             &  source_up_clear, source_dn_clear)
257        ! Use adding method to compute fluxes
258        call adding_ica_lw(ng, nlev, &
259             &  ref_clear, trans_clear, source_up_clear, source_dn_clear, &
260             &  emission(:,jcol), albedo(:,jcol), &
261             &  flux_up_clear, flux_dn_clear)
262      end if
263
264      if (config%do_clear) then
265        ! Sum over g-points to compute broadband fluxes
266        do jlev = 1,nlev+1
267          sum_up = 0.0_jprb
268          sum_dn = 0.0_jprb
269          !$omp simd reduction(+:sum_up, sum_dn)
270          do jg = 1,ng
271            sum_up = sum_up + flux_up_clear(jg,jlev)
272            sum_dn = sum_dn + flux_dn_clear(jg,jlev)
273          end do
274          flux%lw_up_clear(jcol,jlev) = sum_up
275          flux%lw_dn_clear(jcol,jlev) = sum_dn
276        end do
277
278        ! Store surface spectral downwelling fluxes / TOA upwelling
279        do jg = 1,ng
280          flux%lw_dn_surf_clear_g(jg,jcol) = flux_dn_clear(jg,nlev+1)
281          flux%lw_up_toa_clear_g (jg,jcol) = flux_up_clear(jg,1)
282        end do
283        ! Save the spectral fluxes if required
284        if (config%do_save_spectral_flux) then
285          do jlev = 1,nlev+1
286            call indexed_sum(flux_up_clear(:,jlev), &
287                 &           config%i_spec_from_reordered_g_lw, &
288                 &           flux%lw_up_clear_band(:,jcol,jlev))
289            call indexed_sum(flux_dn_clear(:,jlev), &
290                 &           config%i_spec_from_reordered_g_lw, &
291                 &           flux%lw_dn_clear_band(:,jcol,jlev))
292          end do
293        end if
294      end if
295
296      ! --------------------------------------------------------
297      ! Section 4: Loop over cloudy layers to compute reflectance and transmittance
298      ! --------------------------------------------------------
299      ! In this section the reflectance, transmittance and sources
300      ! are computed for each layer
301     
302      ! Firstly, ensure clear-sky transmittance is valid for whole
303      ! depth of the atmosphere, because even above cloud it is used
304      ! by the LW derivatives
305      transmittance(:,1,:) = trans_clear(:,:)
306      ! Dummy values in cloudy regions above cloud top
307      if (i_cloud_top > 0) then
308        transmittance(:,2:,1:min(i_cloud_top,nlev)) = 1.0_jprb
309      end if
310
311      do jlev = i_cloud_top,nlev ! Start at cloud top and work down
312
313        ! Copy over clear-sky properties
314        reflectance(:,1,jlev)    = ref_clear(:,jlev)
315        source_up(:,1,jlev)      = source_up_clear(:,jlev) ! Scaled later by region size
316        source_dn(:,1,jlev)      = source_dn_clear(:,jlev) ! Scaled later by region size
317        nreg = nregions
318        if (is_clear_sky_layer(jlev)) then
319          nreg = 1
320          reflectance(:,2:,jlev)   = 0.0_jprb
321          transmittance(:,2:,jlev) = 1.0_jprb
322          source_up(:,2:,jlev)     = 0.0_jprb
323          source_dn(:,2:,jlev)     = 0.0_jprb
324        else
325          do jreg = 2,nreg
326            ! Cloudy sky
327            ! Add scaled cloud optical depth to clear-sky value
328            od_cloud_new = od_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
329                 &       * od_scaling(jreg,jlev,jcol)
330            od_total = od(:,jlev,jcol) + od_cloud_new
331
332            if (config%do_lw_cloud_scattering) then
333              ssa_total = 0.0_jprb
334              g_total   = 0.0_jprb             
335              if (config%do_lw_aerosol_scattering) then
336                where (od_total > 0.0_jprb)
337                  ssa_total = (ssa(:,jlev,jcol)*od(:,jlev,jcol) &
338                       &     + ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
339                       &     *  od_cloud_new) &
340                       &     / od_total
341                end where
342                where (ssa_total > 0.0_jprb .and. od_total > 0.0_jprb)
343                  g_total = (g(:,jlev,jcol)*ssa(:,jlev,jcol)*od(:,jlev,jcol) &
344                       &     +   g_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
345                       &     * ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
346                       &     *  od_cloud_new) &
347                       &     / (ssa_total*od_total)
348                end where
349              else
350                where (od_total > 0.0_jprb)
351                  ssa_total = ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
352                       &     * od_cloud_new / od_total
353                end where
354                where (ssa_total > 0.0_jprb .and. od_total > 0.0_jprb)
355                  g_total = g_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
356                       &     * ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
357                       &     *  od_cloud_new / (ssa_total*od_total)
358                end where
359              end if
360              call calc_ref_trans_lw(ng, &
361                   &  od_total, ssa_total, g_total, &
362                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1,jcol), &
363                   &  reflectance(:,jreg,jlev), transmittance(:,jreg,jlev), &
364                   &  source_up(:,jreg,jlev), source_dn(:,jreg,jlev))
365            else
366              ! No-scattering case: use simpler functions for
367              ! transmission and emission
368              call calc_no_scattering_transmittance_lw(ng, od_total, &
369                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1, jcol), &
370                   &  transmittance(:,jreg,jlev), source_up(:,jreg,jlev), source_dn(:,jreg,jlev))
371              reflectance(:,jreg,jlev) = 0.0_jprb
372            end if
373          end do
374          ! Emission is scaled by the size of each region
375          do jreg = 1,nregions
376            source_up(:,jreg,jlev) = region_fracs(jreg,jlev,jcol) * source_up(:,jreg,jlev)
377            source_dn(:,jreg,jlev) = region_fracs(jreg,jlev,jcol) * source_dn(:,jreg,jlev)
378          end do
379        end if
380
381      end do ! Loop over levels
382
383      ! --------------------------------------------------------
384      ! Section 5: Compute total sources and albedos at each half level
385      ! --------------------------------------------------------
386
387      total_albedo(:,:,:) = 0.0_jprb
388      total_source(:,:,:) = 0.0_jprb
389
390      ! Calculate the upwelling radiation emitted by the surface, and
391      ! copy the surface albedo into total_albedo
392      do jreg = 1,nregions
393        do jg = 1,ng
394          ! region_fracs(jreg,nlev,jcol) is the fraction of each region in the
395          ! lowest model level
396          total_source(jg,jreg,nlev+1) = region_fracs(jreg,nlev,jcol)*emission(jg,jcol)
397          total_albedo(jg,jreg,nlev+1) = albedo(jg,jcol)
398        end do
399      end do
400
401      ! Work up from the surface computing the total albedo of the
402      ! atmosphere and the total upwelling due to emission below each
403      ! level below using the adding method
404      do jlev = nlev,i_cloud_top,-1
405
406        total_albedo_below        = 0.0_jprb
407
408        if (is_clear_sky_layer(jlev)) then
409          total_albedo_below = 0.0_jprb
410          total_source_below = 0.0_jprb
411          do jg = 1,ng
412            inv_denom(jg,1) = 1.0_jprb &
413                 &  / (1.0_jprb - total_albedo(jg,1,jlev+1)*reflectance(jg,1,jlev))
414            total_albedo_below(jg,1) = reflectance(jg,1,jlev) &
415                 &  + transmittance(jg,1,jlev)*transmittance(jg,1,jlev)*total_albedo(jg,1,jlev+1) &
416                 &  * inv_denom(jg,1)
417            total_source_below(jg,1) = source_up(jg,1,jlev) &
418                 &  + transmittance(jg,1,jlev)*(total_source(jg,1,jlev+1) &
419                 &  + total_albedo(jg,1,jlev+1)*source_dn(jg,1,jlev)) &
420                 &  * inv_denom(jg,1)
421          end do
422        else
423          inv_denom = 1.0_jprb / (1.0_jprb - total_albedo(:,:,jlev+1)*reflectance(:,:,jlev))
424          total_albedo_below = reflectance(:,:,jlev) &
425               &  + transmittance(:,:,jlev)*transmittance(:,:,jlev)*total_albedo(:,:,jlev+1) &
426               &  * inv_denom
427          total_source_below = source_up(:,:,jlev) &
428               &  + transmittance(:,:,jlev)*(total_source(:,:,jlev+1) &
429               &  + total_albedo(:,:,jlev+1)*source_dn(:,:,jlev)) &
430               &  * inv_denom
431        end if
432
433        ! Account for cloud overlap when converting albedo below a
434        ! layer interface to the equivalent values just above
435        if (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev-1)) then
436          total_albedo(:,:,jlev) = total_albedo_below(:,:)
437          total_source(:,:,jlev) = total_source_below(:,:)
438        else
439          total_source(:,:,jlev) = singlemat_x_vec(ng,ng,nregions,&
440               &  u_matrix(:,:,jlev,jcol), total_source_below)
441          ! Use overlap matrix and exclude "anomalous" horizontal
442          ! transport described by Shonk & Hogan (2008).  Therefore,
443          ! the operation we perform is essentially diag(total_albedo)
444          ! = matmul(transpose(v_matrix), diag(total_albedo_below)).
445          do jreg = 1,nregions
446            do jreg2 = 1,nregions
447              total_albedo(:,jreg,jlev) &
448                   &  = total_albedo(:,jreg,jlev) &
449                   &  + total_albedo_below(:,jreg2) &
450                   &  * v_matrix(jreg2,jreg,jlev,jcol)
451
452            end do
453          end do
454
455        end if
456       
457      end do ! Reverse loop over levels
458
459      ! --------------------------------------------------------
460      ! Section 6: Copy over downwelling fluxes above cloud top
461      ! --------------------------------------------------------
462      do jlev = 1,i_cloud_top
463        if (config%do_clear) then
464          ! Clear-sky fluxes have already been averaged: use these
465          flux%lw_dn(jcol,jlev) = flux%lw_dn_clear(jcol,jlev)
466          if (config%do_save_spectral_flux) then
467            flux%lw_dn_band(:,jcol,jlev) = flux%lw_dn_clear_band(:,jcol,jlev)
468          end if
469        else
470          sum_dn = 0.0_jprb
471          !$omp simd reduction(+:sum_dn)
472          do jg = 1,ng
473            sum_dn = sum_dn + flux_dn_clear(jg,jlev)
474          end do
475          flux%lw_dn(jcol,jlev) = sum_dn
476          if (config%do_save_spectral_flux) then
477            call indexed_sum(flux_dn_clear(:,jlev), &
478                 &           config%i_spec_from_reordered_g_lw, &
479                 &           flux%lw_dn_band(:,jcol,jlev))
480          end if
481        end if
482      end do
483
484      ! --------------------------------------------------------
485      ! Section 7: Compute fluxes up to top-of-atmosphere
486      ! --------------------------------------------------------
487
488      ! Compute the fluxes just above the highest cloud
489      flux_up(:,1) = total_source(:,1,i_cloud_top) &
490           &  + total_albedo(:,1,i_cloud_top)*flux_dn_clear(:,i_cloud_top)
491      flux_up(:,2:) = 0.0_jprb
492
493      sum_up = 0.0_jprb
494      !$omp simd reduction(+:sum_up)
495      do jg = 1,ng
496        sum_up = sum_up + flux_up(jg,1)
497      end do
498      flux%lw_up(jcol,i_cloud_top) = sum_up
499
500      if (config%do_save_spectral_flux) then
501        call indexed_sum(flux_up(:,1), &
502             &           config%i_spec_from_reordered_g_lw, &
503             &           flux%lw_up_band(:,jcol,i_cloud_top))
504      end if
505      do jlev = i_cloud_top-1,1,-1
506        flux_up(:,1) = trans_clear(:,jlev)*flux_up(:,1) + source_up_clear(:,jlev)
507        sum_up = 0.0_jprb
508        !$omp simd reduction(+:sum_up)
509        do jg = 1,ng
510          sum_up = sum_up + flux_up(jg,1)
511        end do
512        flux%lw_up(jcol,jlev) = sum_up
513        if (config%do_save_spectral_flux) then
514          call indexed_sum(flux_up(:,1), &
515               &           config%i_spec_from_reordered_g_lw, &
516               &           flux%lw_up_band(:,jcol,jlev))
517        end if
518      end do
519      flux%lw_up_toa_g(:,jcol) = sum(flux_up,2)
520
521      ! --------------------------------------------------------
522      ! Section 8: Compute fluxes down to surface
523      ! --------------------------------------------------------
524
525      ! Copy over downwelling spectral fluxes at top of first
526      ! scattering layer, using overlap matrix to translate to the
527      ! regions of the first layer of cloud
528      do jreg = 1,nregions
529        flux_dn(:,jreg)  = v_matrix(jreg,1,i_cloud_top,jcol)*flux_dn_clear(:,i_cloud_top)
530      end do
531
532      ! Final loop back down through the atmosphere to compute fluxes
533      do jlev = i_cloud_top,nlev
534
535        if (is_clear_sky_layer(jlev)) then
536          do jg = 1,ng
537            flux_dn(jg,1) = (transmittance(jg,1,jlev)*flux_dn(jg,1) &
538                 &  + reflectance(jg,1,jlev)*total_source(jg,1,jlev+1) + source_dn(jg,1,jlev) ) &
539                 &  / (1.0_jprb - reflectance(jg,1,jlev)*total_albedo(jg,1,jlev+1))
540            flux_up(jg,1) = total_source(jg,1,jlev+1) + flux_dn(jg,1)*total_albedo(jg,1,jlev+1)
541          end do
542          flux_dn(:,2:)  = 0.0_jprb
543          flux_up(:,2:)  = 0.0_jprb
544        else
545          flux_dn = (transmittance(:,:,jlev)*flux_dn &
546               &     + reflectance(:,:,jlev)*total_source(:,:,jlev+1) + source_dn(:,:,jlev) ) &
547               &  / (1.0_jprb - reflectance(:,:,jlev)*total_albedo(:,:,jlev+1))
548          flux_up = total_source(:,:,jlev+1) + flux_dn*total_albedo(:,:,jlev+1)
549        end if
550
551        if (.not. (is_clear_sky_layer(jlev) &
552             &    .and. is_clear_sky_layer(jlev+1))) then
553          ! Account for overlap rules in translating fluxes just above
554          ! a layer interface to the values just below
555          flux_dn_below = singlemat_x_vec(ng,ng,nregions, &
556               &  v_matrix(:,:,jlev+1,jcol), flux_dn)
557          flux_dn = flux_dn_below
558        end if ! Otherwise the fluxes in each region are the same so
559               ! nothing to do
560
561        ! Store the broadband fluxes
562        sum_up = 0.0_jprb
563        sum_dn = 0.0_jprb
564        do jreg = 1,nregions
565          !$omp simd reduction(+:sum_up, sum_dn)
566          do jg = 1,ng
567            sum_up = sum_up + flux_up(jg,jreg)
568            sum_dn = sum_dn + flux_dn(jg,jreg)
569          end do
570        end do
571        flux%lw_up(jcol,jlev+1) = sum_up
572        flux%lw_dn(jcol,jlev+1) = sum_dn
573
574        ! Save the spectral fluxes if required
575        if (config%do_save_spectral_flux) then
576          call indexed_sum(sum(flux_up,2), &
577               &           config%i_spec_from_reordered_g_lw, &
578               &           flux%lw_up_band(:,jcol,jlev+1))
579          call indexed_sum(sum(flux_dn,2), &
580               &           config%i_spec_from_reordered_g_lw, &
581               &           flux%lw_dn_band(:,jcol,jlev+1))
582         end if
583
584      end do ! Final loop over levels
585     
586      ! Store surface spectral downwelling fluxes, which at this point
587      ! are at the surface
588      flux%lw_dn_surf_g(:,jcol) = sum(flux_dn,2)
589
590      ! Compute the longwave derivatives needed by Hogan and Bozzo
591      ! (2015) approximate radiation update scheme
592      if (config%do_lw_derivatives) then
593        ! Note that at this point flux_up contains the spectral
594        ! fluxes into the regions of the lowest layer; we sum over
595        ! regions first to provide a simple spectral flux upwelling
596        ! from the surface
597        call calc_lw_derivatives_region(ng, nlev, nregions, jcol, transmittance, &
598             &  u_matrix(:,:,:,jcol), sum(flux_up,2), flux%lw_derivatives)
599      end if
600     
601    end do ! Loop over columns
602
603    if (lhook) call dr_hook('radiation_tripleclouds_lw:solver_tripleclouds_lw',1,hook_handle)
604
605  end subroutine solver_tripleclouds_lw
606
607end module radiation_tripleclouds_lw
Note: See TracBrowser for help on using the repository browser.