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

Last change on this file since 4773 was 4773, checked in by idelkadi, 5 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


File size: 24.5 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    ! Index of the highest cloudy layer
173    integer :: i_cloud_top
174
175    integer :: jcol, jlev, jg, jreg, jreg2, ng
176
177    real(jphook) :: hook_handle
178
179    if (lhook) call dr_hook('radiation_tripleclouds_lw:solver_tripleclouds_lw',0,hook_handle)
180
181    ! --------------------------------------------------------
182    ! Section 1: Prepare general variables and arrays
183    ! --------------------------------------------------------
184    ! Copy array dimensions to local variables for convenience
185    ng   = config%n_g_lw
186
187    ! Compute the wavelength-independent region fractions and
188    ! optical-depth scalings
189    call calc_region_properties(nlev,nregions,istartcol,iendcol, &
190         &  config%i_cloud_pdf_shape == IPdfShapeGamma, &
191         &  cloud%fraction, cloud%fractional_std, region_fracs, &
192         &  od_scaling, config%cloud_fraction_threshold)
193
194    ! Compute wavelength-independent overlap matrices u_matrix and v_matrix
195    call calc_overlap_matrices(nlev,nregions,istartcol,iendcol, &
196         &  region_fracs, cloud%overlap_param, &
197         &  u_matrix, v_matrix, &
198         &  decorrelation_scaling=config%cloud_inhom_decorr_scaling, &
199         &  cloud_fraction_threshold=config%cloud_fraction_threshold, &
200         &  use_beta_overlap=config%use_beta_overlap, &
201         &  cloud_cover=flux%cloud_cover_lw)
202
203    ! Main loop over columns
204    do jcol = istartcol, iendcol
205      ! --------------------------------------------------------
206      ! Section 2: Prepare column-specific variables and arrays
207      ! --------------------------------------------------------
208
209      ! Define which layers contain cloud; assume that
210      ! cloud%crop_cloud_fraction has already been called
211      is_clear_sky_layer = .true.
212      i_cloud_top = nlev+1
213      do jlev = 1,nlev
214        if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
215          is_clear_sky_layer(jlev) = .false.
216          ! Get index to the first cloudy layer from the top
217          if (i_cloud_top > jlev) then
218            i_cloud_top = jlev
219          end if
220        end if
221      end do
222      if (config%do_lw_aerosol_scattering) then
223        ! This is actually the first layer in which we need to
224        ! consider scattering
225        i_cloud_top = 1
226      end if
227
228      ! --------------------------------------------------------
229      ! Section 3: Clear-sky calculation
230      ! --------------------------------------------------------
231
232      if (.not. config%do_lw_aerosol_scattering) then
233        ! No scattering in clear-sky flux calculation; note that here
234        ! the first two dimensions of the input arrays are unpacked
235        ! into vectors inside the routine       
236        call calc_no_scattering_transmittance_lw(ng*nlev, od(:,:,jcol), &
237             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1, jcol), &
238             &  trans_clear, source_up_clear, source_dn_clear)
239        ! Ensure that clear-sky reflectance is zero since it may be
240        ! used in cloudy-sky case
241        ref_clear = 0.0_jprb
242        ! Simple down-then-up method to compute fluxes
243        call calc_fluxes_no_scattering_lw(ng, nlev, &
244             &  trans_clear, source_up_clear, source_dn_clear, &
245             &  emission(:,jcol), albedo(:,jcol), &
246             &  flux_up_clear, flux_dn_clear)
247      else
248        ! Scattering in clear-sky flux calculation
249        call calc_ref_trans_lw(ng*nlev, &
250             &  od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), &
251             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1,jcol), &
252             &  ref_clear, trans_clear, &
253             &  source_up_clear, source_dn_clear)
254        ! Use adding method to compute fluxes
255        call adding_ica_lw(ng, nlev, &
256             &  ref_clear, trans_clear, source_up_clear, source_dn_clear, &
257             &  emission(:,jcol), albedo(:,jcol), &
258             &  flux_up_clear, flux_dn_clear)
259      end if
260
261      if (config%do_clear) then
262        ! Sum over g-points to compute broadband fluxes
263        flux%lw_up_clear(jcol,:) = sum(flux_up_clear,1)
264        flux%lw_dn_clear(jcol,:) = sum(flux_dn_clear,1)
265        ! Store surface spectral downwelling fluxes / TOA upwelling
266        flux%lw_dn_surf_clear_g(:,jcol) = flux_dn_clear(:,nlev+1)
267        flux%lw_up_toa_clear_g (:,jcol) = flux_up_clear(:,1)
268        ! Save the spectral fluxes if required
269        if (config%do_save_spectral_flux) then
270          do jlev = 1,nlev+1
271            call indexed_sum(flux_up_clear(:,jlev), &
272                 &           config%i_spec_from_reordered_g_lw, &
273                 &           flux%lw_up_clear_band(:,jcol,jlev))
274            call indexed_sum(flux_dn_clear(:,jlev), &
275                 &           config%i_spec_from_reordered_g_lw, &
276                 &           flux%lw_dn_clear_band(:,jcol,jlev))
277          end do
278        end if
279      end if
280
281      ! --------------------------------------------------------
282      ! Section 4: Loop over cloudy layers to compute reflectance and transmittance
283      ! --------------------------------------------------------
284      ! In this section the reflectance, transmittance and sources
285      ! are computed for each layer
286     
287      ! Firstly, ensure clear-sky transmittance is valid for whole
288      ! depth of the atmosphere, because even above cloud it is used
289      ! by the LW derivatives
290      transmittance(:,1,:) = trans_clear(:,:)
291      ! Dummy values in cloudy regions above cloud top
292      if (i_cloud_top > 0) then
293        transmittance(:,2:,1:min(i_cloud_top,nlev)) = 1.0_jprb
294      end if
295
296      do jlev = i_cloud_top,nlev ! Start at cloud top and work down
297
298        ! Copy over clear-sky properties
299        reflectance(:,1,jlev)    = ref_clear(:,jlev)
300        source_up(:,1,jlev)      = source_up_clear(:,jlev) ! Scaled later by region size
301        source_dn(:,1,jlev)      = source_dn_clear(:,jlev) ! Scaled later by region size
302        nreg = nregions
303        if (is_clear_sky_layer(jlev)) then
304          nreg = 1
305          reflectance(:,2:,jlev)   = 0.0_jprb
306          transmittance(:,2:,jlev) = 1.0_jprb
307          source_up(:,2:,jlev)     = 0.0_jprb
308          source_dn(:,2:,jlev)     = 0.0_jprb
309        else
310          do jreg = 2,nreg
311            ! Cloudy sky
312            ! Add scaled cloud optical depth to clear-sky value
313            od_cloud_new = od_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
314                 &       * od_scaling(jreg,jlev,jcol)
315            od_total = od(:,jlev,jcol) + od_cloud_new
316
317            if (config%do_lw_cloud_scattering) then
318              ssa_total = 0.0_jprb
319              g_total   = 0.0_jprb             
320              if (config%do_lw_aerosol_scattering) then
321                where (od_total > 0.0_jprb)
322                  ssa_total = (ssa(:,jlev,jcol)*od(:,jlev,jcol) &
323                       &     + ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
324                       &     *  od_cloud_new) &
325                       &     / od_total
326                end where
327                where (ssa_total > 0.0_jprb .and. od_total > 0.0_jprb)
328                  g_total = (g(:,jlev,jcol)*ssa(:,jlev,jcol)*od(:,jlev,jcol) &
329                       &     +   g_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
330                       &     * ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
331                       &     *  od_cloud_new) &
332                       &     / (ssa_total*od_total)
333                end where
334              else
335                where (od_total > 0.0_jprb)
336                  ssa_total = ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
337                       &     * od_cloud_new / od_total
338                end where
339                where (ssa_total > 0.0_jprb .and. od_total > 0.0_jprb)
340                  g_total = g_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
341                       &     * ssa_cloud(config%i_band_from_reordered_g_lw,jlev,jcol) &
342                       &     *  od_cloud_new / (ssa_total*od_total)
343                end where
344              end if
345              call calc_ref_trans_lw(ng, &
346                   &  od_total, ssa_total, g_total, &
347                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1,jcol), &
348                   &  reflectance(:,jreg,jlev), transmittance(:,jreg,jlev), &
349                   &  source_up(:,jreg,jlev), source_dn(:,jreg,jlev))
350            else
351              ! No-scattering case: use simpler functions for
352              ! transmission and emission
353              call calc_no_scattering_transmittance_lw(ng, od_total, &
354                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1, jcol), &
355                   &  transmittance(:,jreg,jlev), source_up(:,jreg,jlev), source_dn(:,jreg,jlev))
356              reflectance(:,jreg,jlev) = 0.0_jprb
357            end if
358          end do
359          ! Emission is scaled by the size of each region
360          do jreg = 1,nregions
361            source_up(:,jreg,jlev) = region_fracs(jreg,jlev,jcol) * source_up(:,jreg,jlev)
362            source_dn(:,jreg,jlev) = region_fracs(jreg,jlev,jcol) * source_dn(:,jreg,jlev)
363          end do
364        end if
365
366      end do ! Loop over levels
367
368      ! --------------------------------------------------------
369      ! Section 5: Compute total sources and albedos at each half level
370      ! --------------------------------------------------------
371
372      total_albedo(:,:,:) = 0.0_jprb
373      total_source(:,:,:) = 0.0_jprb
374
375      ! Calculate the upwelling radiation emitted by the surface, and
376      ! copy the surface albedo into total_albedo
377      do jreg = 1,nregions
378        do jg = 1,ng
379          ! region_fracs(jreg,nlev,jcol) is the fraction of each region in the
380          ! lowest model level
381          total_source(jg,jreg,nlev+1) = region_fracs(jreg,nlev,jcol)*emission(jg,jcol)
382          total_albedo(jg,jreg,nlev+1) = albedo(jg,jcol)
383        end do
384      end do
385
386      ! Work up from the surface computing the total albedo of the
387      ! atmosphere and the total upwelling due to emission below each
388      ! level below using the adding method
389      do jlev = nlev,i_cloud_top,-1
390
391        total_albedo_below        = 0.0_jprb
392
393        if (is_clear_sky_layer(jlev)) then
394          total_albedo_below = 0.0_jprb
395          total_source_below = 0.0_jprb
396          do jg = 1,ng
397            inv_denom(jg,1) = 1.0_jprb &
398                 &  / (1.0_jprb - total_albedo(jg,1,jlev+1)*reflectance(jg,1,jlev))
399            total_albedo_below(jg,1) = reflectance(jg,1,jlev) &
400                 &  + transmittance(jg,1,jlev)*transmittance(jg,1,jlev)*total_albedo(jg,1,jlev+1) &
401                 &  * inv_denom(jg,1)
402            total_source_below(jg,1) = source_up(jg,1,jlev) &
403                 &  + transmittance(jg,1,jlev)*(total_source(jg,1,jlev+1) &
404                 &  + total_albedo(jg,1,jlev+1)*source_dn(jg,1,jlev)) &
405                 &  * inv_denom(jg,1)
406          end do
407        else
408          inv_denom = 1.0_jprb / (1.0_jprb - total_albedo(:,:,jlev+1)*reflectance(:,:,jlev))
409          total_albedo_below = reflectance(:,:,jlev) &
410               &  + transmittance(:,:,jlev)*transmittance(:,:,jlev)*total_albedo(:,:,jlev+1) &
411               &  * inv_denom
412          total_source_below = source_up(:,:,jlev) &
413               &  + transmittance(:,:,jlev)*(total_source(:,:,jlev+1) &
414               &  + total_albedo(:,:,jlev+1)*source_dn(:,:,jlev)) &
415               &  * inv_denom
416        end if
417
418        ! Account for cloud overlap when converting albedo below a
419        ! layer interface to the equivalent values just above
420        if (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev-1)) then
421          total_albedo(:,:,jlev) = total_albedo_below(:,:)
422          total_source(:,:,jlev) = total_source_below(:,:)
423        else
424          total_source(:,:,jlev) = singlemat_x_vec(ng,ng,nregions,&
425               &  u_matrix(:,:,jlev,jcol), total_source_below)
426          ! Use overlap matrix and exclude "anomalous" horizontal
427          ! transport described by Shonk & Hogan (2008).  Therefore,
428          ! the operation we perform is essentially diag(total_albedo)
429          ! = matmul(transpose(v_matrix), diag(total_albedo_below)).
430          do jreg = 1,nregions
431            do jreg2 = 1,nregions
432              total_albedo(:,jreg,jlev) &
433                   &  = total_albedo(:,jreg,jlev) &
434                   &  + total_albedo_below(:,jreg2) &
435                   &  * v_matrix(jreg2,jreg,jlev,jcol)
436
437            end do
438          end do
439
440        end if
441       
442      end do ! Reverse loop over levels
443
444      ! --------------------------------------------------------
445      ! Section 6: Copy over downwelling fluxes above cloud top
446      ! --------------------------------------------------------
447      do jlev = 1,i_cloud_top
448        if (config%do_clear) then
449          ! Clear-sky fluxes have already been averaged: use these
450          flux%lw_dn(jcol,jlev) = flux%lw_dn_clear(jcol,jlev)
451          if (config%do_save_spectral_flux) then
452            flux%lw_dn_band(:,jcol,jlev) = flux%lw_dn_clear_band(:,jcol,jlev)
453          end if
454        else
455          flux%lw_dn(jcol,:) = sum(flux_dn_clear(:,jlev))
456          if (config%do_save_spectral_flux) then
457            call indexed_sum(flux_dn_clear(:,jlev), &
458                 &           config%i_spec_from_reordered_g_lw, &
459                 &           flux%lw_dn_band(:,jcol,jlev))
460          end if
461        end if
462      end do
463
464      ! --------------------------------------------------------
465      ! Section 7: Compute fluxes up to top-of-atmosphere
466      ! --------------------------------------------------------
467
468      ! Compute the fluxes just above the highest cloud
469      flux_up(:,1) = total_source(:,1,i_cloud_top) &
470           &  + total_albedo(:,1,i_cloud_top)*flux_dn_clear(:,i_cloud_top)
471      flux_up(:,2:) = 0.0_jprb
472      flux%lw_up(jcol,i_cloud_top) = sum(flux_up(:,1))
473      if (config%do_save_spectral_flux) then
474        call indexed_sum(flux_up(:,1), &
475             &           config%i_spec_from_reordered_g_lw, &
476             &           flux%lw_up_band(:,jcol,i_cloud_top))
477      end if
478      do jlev = i_cloud_top-1,1,-1
479        flux_up(:,1) = trans_clear(:,jlev)*flux_up(:,1) + source_up_clear(:,jlev)
480        flux%lw_up(jcol,jlev) = sum(flux_up(:,1))
481        if (config%do_save_spectral_flux) then
482          call indexed_sum(flux_up(:,1), &
483               &           config%i_spec_from_reordered_g_lw, &
484               &           flux%lw_up_band(:,jcol,jlev))
485        end if
486      end do
487      flux%lw_up_toa_g(:,jcol) = sum(flux_up,2)
488
489      ! --------------------------------------------------------
490      ! Section 8: Compute fluxes down to surface
491      ! --------------------------------------------------------
492
493      ! Copy over downwelling spectral fluxes at top of first
494      ! scattering layer, using overlap matrix to translate to the
495      ! regions of the first layer of cloud
496      do jreg = 1,nregions
497        flux_dn(:,jreg)  = v_matrix(jreg,1,i_cloud_top,jcol)*flux_dn_clear(:,i_cloud_top)
498      end do
499
500      ! Final loop back down through the atmosphere to compute fluxes
501      do jlev = i_cloud_top,nlev
502
503        if (is_clear_sky_layer(jlev)) then
504          do jg = 1,ng
505            flux_dn(jg,1) = (transmittance(jg,1,jlev)*flux_dn(jg,1) &
506                 &  + reflectance(jg,1,jlev)*total_source(jg,1,jlev+1) + source_dn(jg,1,jlev) ) &
507                 &  / (1.0_jprb - reflectance(jg,1,jlev)*total_albedo(jg,1,jlev+1))
508            flux_up(jg,1) = total_source(jg,1,jlev+1) + flux_dn(jg,1)*total_albedo(jg,1,jlev+1)
509          end do
510          flux_dn(:,2:)  = 0.0_jprb
511          flux_up(:,2:)  = 0.0_jprb
512        else
513          flux_dn = (transmittance(:,:,jlev)*flux_dn &
514               &     + reflectance(:,:,jlev)*total_source(:,:,jlev+1) + source_dn(:,:,jlev) ) &
515               &  / (1.0_jprb - reflectance(:,:,jlev)*total_albedo(:,:,jlev+1))
516          flux_up = total_source(:,:,jlev+1) + flux_dn*total_albedo(:,:,jlev+1)
517        end if
518
519        if (.not. (is_clear_sky_layer(jlev) &
520             &    .and. is_clear_sky_layer(jlev+1))) then
521          ! Account for overlap rules in translating fluxes just above
522          ! a layer interface to the values just below
523          flux_dn_below = singlemat_x_vec(ng,ng,nregions, &
524               &  v_matrix(:,:,jlev+1,jcol), flux_dn)
525          flux_dn = flux_dn_below
526        end if ! Otherwise the fluxes in each region are the same so
527               ! nothing to do
528
529        ! Store the broadband fluxes
530        flux%lw_up(jcol,jlev+1) = sum(sum(flux_up,1))
531        flux%lw_dn(jcol,jlev+1) = sum(sum(flux_dn,1))
532
533        ! Save the spectral fluxes if required
534        if (config%do_save_spectral_flux) then
535          call indexed_sum(sum(flux_up,2), &
536               &           config%i_spec_from_reordered_g_lw, &
537               &           flux%lw_up_band(:,jcol,jlev+1))
538          call indexed_sum(sum(flux_dn,2), &
539               &           config%i_spec_from_reordered_g_lw, &
540               &           flux%lw_dn_band(:,jcol,jlev+1))
541         end if
542
543      end do ! Final loop over levels
544     
545      ! Store surface spectral downwelling fluxes, which at this point
546      ! are at the surface
547      flux%lw_dn_surf_g(:,jcol) = sum(flux_dn,2)
548
549      ! Compute the longwave derivatives needed by Hogan and Bozzo
550      ! (2015) approximate radiation update scheme
551      if (config%do_lw_derivatives) then
552        ! Note that at this point flux_up contains the spectral
553        ! fluxes into the regions of the lowest layer; we sum over
554        ! regions first to provide a simple spectral flux upwelling
555        ! from the surface
556        call calc_lw_derivatives_region(ng, nlev, nregions, jcol, transmittance, &
557             &  u_matrix(:,:,:,jcol), sum(flux_up,2), flux%lw_derivatives)
558      end if
559     
560    end do ! Loop over columns
561
562    if (lhook) call dr_hook('radiation_tripleclouds_lw:solver_tripleclouds_lw',1,hook_handle)
563
564  end subroutine solver_tripleclouds_lw
565
566end module radiation_tripleclouds_lw
Note: See TracBrowser for help on using the repository browser.