source: LMDZ6/trunk/libf/phylmd/ecrad/radiation_tripleclouds_sw.F90 @ 3908

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

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 25.9 KB
Line 
1! radiation_tripleclouds_sw.F90 - Shortwave "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-11  R. Hogan  Receive albedos at g-points
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!   2019-01-02  R. Hogan  Fixed problem of do_save_spectral_flux .and. .not. do_sw_direct
21
22module radiation_tripleclouds_sw
23
24  public
25
26contains
27  ! Provides elemental function "delta_eddington"
28#include "radiation_delta_eddington.h"
29
30  ! Small routine for scaling cloud optical depth in the cloudy
31  ! regions
32#include "radiation_optical_depth_scaling.h"
33
34  ! This module contains just one subroutine, the shortwave
35  ! "Tripleclouds" solver in which cloud inhomogeneity is treated by
36  ! dividing each model level into three regions, one clear and two
37  ! cloudy (with differing optical depth). This approach was described
38  ! by Shonk and Hogan (2008).
39
40  subroutine solver_tripleclouds_sw(nlev,istartcol,iendcol, &
41       &  config, single_level, cloud, &
42       &  od, ssa, g, od_cloud, ssa_cloud, g_cloud, &
43       &  albedo_direct, albedo_diffuse, incoming_sw, &
44       &  flux)
45
46    use parkind1, only           : jprb
47    use yomhook,  only           : lhook, dr_hook
48
49!    use radiation_io, only             : nulout
50    use radiation_config, only         : config_type, IPdfShapeGamma
51    use radiation_single_level, only   : single_level_type
52    use radiation_cloud, only          : cloud_type
53    use radiation_regions, only        : calc_region_properties
54    use radiation_overlap, only        : calc_overlap_matrices
55    use radiation_flux, only           : flux_type, &
56         &                               indexed_sum, add_indexed_sum
57    use radiation_matrix, only         : singlemat_x_vec
58    use radiation_two_stream, only     : calc_two_stream_gammas_sw, &
59         &                       calc_reflectance_transmittance_sw
60
61    implicit none
62
63    ! Inputs
64    integer, intent(in) :: nlev               ! number of model levels
65    integer, intent(in) :: istartcol, iendcol ! range of columns to process
66    type(config_type),        intent(in) :: config
67    type(single_level_type),  intent(in) :: single_level
68    type(cloud_type),         intent(in) :: cloud
69
70    ! Gas and aerosol optical depth, single-scattering albedo and
71    ! asymmetry factor at each shortwave g-point
72!    real(jprb), intent(in), dimension(istartcol:iendcol,nlev,config%n_g_sw) :: &
73    real(jprb), intent(in), dimension(config%n_g_sw,nlev,istartcol:iendcol) :: &
74         &  od, ssa, g
75
76    ! Cloud and precipitation optical depth, single-scattering albedo and
77    ! asymmetry factor in each shortwave band
78    real(jprb), intent(in), dimension(config%n_bands_sw,nlev,istartcol:iendcol) :: &
79         &  od_cloud, ssa_cloud, g_cloud
80
81    ! Optical depth, single scattering albedo and asymmetry factor in
82    ! each g-point including gas, aerosol and clouds
83    real(jprb), dimension(config%n_g_sw) :: od_total, ssa_total, g_total
84
85    ! Direct and diffuse surface albedos, and the incoming shortwave
86    ! flux into a plane perpendicular to the incoming radiation at
87    ! top-of-atmosphere in each of the shortwave g points
88    real(jprb), intent(in), dimension(config%n_g_sw,istartcol:iendcol) :: &
89         &  albedo_direct, albedo_diffuse, incoming_sw
90
91    ! Output
92    type(flux_type), intent(inout):: flux
93
94    ! Local constants
95    integer, parameter :: nregions = 3
96
97    ! In a clear-sky layer this will be 1, otherwise equal to nregions
98    integer :: nreg
99
100    ! Local variables
101   
102    ! The area fractions of each region
103    real(jprb) :: region_fracs(1:nregions,nlev,istartcol:iendcol)
104
105    ! The scaling used for the optical depth in the cloudy regions
106    real(jprb) :: od_scaling(2:nregions,nlev,istartcol:iendcol)
107
108    ! Directional overlap matrices defined at all layer interfaces
109    ! including top-of-atmosphere and the surface
110    real(jprb), dimension(nregions,nregions,nlev+1, &
111         &                istartcol:iendcol) :: u_matrix, v_matrix
112
113    ! Two-stream variables
114    real(jprb), dimension(config%n_g_sw) :: gamma1, gamma2, gamma3
115
116    ! Diffuse reflection and transmission matrices of each layer
117    real(jprb), dimension(config%n_g_sw, nregions, nlev) &
118         &  :: reflectance, transmittance
119
120    ! Terms translating the direct flux entering the layer from above
121    ! to the reflected radiation exiting upwards (ref_dir) and the
122    ! scattered radiation exiting downwards (trans_dir_diff), along with the
123    ! direct unscattered transmission matrix (trans_dir_dir).
124    real(jprb), dimension(config%n_g_sw, nregions, nlev) &
125         &  :: ref_dir, trans_dir_diff, trans_dir_dir
126
127    ! Total albedo of the atmosphere/surface just above a layer
128    ! interface with respect to downwelling diffuse and direct
129    ! (respectively) radiation at that interface, where level index =
130    ! 1 corresponds to the top-of-atmosphere
131    real(jprb), dimension(config%n_g_sw, nregions, nlev+1) &
132         &  :: total_albedo, total_albedo_direct
133
134    ! ...equivalent values for clear-skies
135    real(jprb), dimension(config%n_g_sw, nlev+1) &
136         &  :: total_albedo_clear, total_albedo_clear_direct
137
138    ! Total albedo of the atmosphere just below a layer interface
139    real(jprb), dimension(config%n_g_sw, nregions) &
140         &  :: total_albedo_below, total_albedo_below_direct
141
142    ! Direct downwelling flux below and above an interface between
143    ! layers into a plane perpendicular to the direction of the sun
144    real(jprb), dimension(config%n_g_sw, nregions) &
145         &  :: direct_dn, direct_dn_below
146    ! Diffuse equivalents
147    real(jprb), dimension(config%n_g_sw, nregions) &
148         &  :: flux_dn, flux_dn_below, flux_up
149
150    ! ...clear-sky equivalent (no distinction between "above/below")
151    real(jprb), dimension(config%n_g_sw) &
152         &  :: direct_dn_clear, flux_dn_clear, flux_up_clear
153
154    ! Clear-sky equivalent, but actually its reciprocal to replace
155    ! some divisions by multiplications
156    real(jprb), dimension(config%n_g_sw, nregions) :: inv_denom
157
158    ! Identify clear-sky layers, with pseudo layers for outer space
159    ! and below the ground, both treated as single-region clear skies
160    logical :: is_clear_sky_layer(0:nlev+1)
161
162    ! Scattering optical depth of gas+aerosol and of cloud
163    real(jprb) :: scat_od, scat_od_cloud
164
165    real(jprb) :: mu0
166
167    integer :: jcol, jlev, jg, jreg, iband, jreg2, ng
168
169    real(jprb) :: hook_handle
170
171    if (lhook) call dr_hook('radiation_tripleclouds_sw:solver_tripleclouds_sw',0,hook_handle)
172
173    ! --------------------------------------------------------
174    ! Section 1: Prepare general variables and arrays
175    ! --------------------------------------------------------
176    ! Copy array dimensions to local variables for convenience
177    ng   = config%n_g_sw
178
179    ! Compute the wavelength-independent region fractions and
180    ! optical-depth scalings
181    call calc_region_properties(nlev,nregions,istartcol,iendcol, &
182         &  config%i_cloud_pdf_shape == IPdfShapeGamma, &
183         &  cloud%fraction, cloud%fractional_std, region_fracs, &
184         &  od_scaling, config%cloud_fraction_threshold)
185
186    ! Compute wavelength-independent overlap matrices u_matrix and
187    ! v_matrix
188    call calc_overlap_matrices(nlev,nregions,istartcol,iendcol, &
189         &  region_fracs, cloud%overlap_param, &
190         &  u_matrix, v_matrix, &
191         &  decorrelation_scaling=config%cloud_inhom_decorr_scaling, &
192         &  cloud_fraction_threshold=config%cloud_fraction_threshold, &
193         &  use_beta_overlap=config%use_beta_overlap, &
194         &  cloud_cover=flux%cloud_cover_sw)
195
196    ! Main loop over columns
197    do jcol = istartcol, iendcol
198      ! --------------------------------------------------------
199      ! Section 2: Prepare column-specific variables and arrays
200      ! --------------------------------------------------------
201
202      ! Copy local cosine of the solar zenith angle
203      mu0 = single_level%cos_sza(jcol)
204
205      ! Skip profile if sun is too low in the sky
206      if (mu0 < 1.0e-10_jprb) then
207        flux%sw_dn(jcol,:) = 0.0_jprb
208        flux%sw_up(jcol,:) = 0.0_jprb
209        if (allocated(flux%sw_dn_direct)) then
210          flux%sw_dn_direct(jcol,:) = 0.0_jprb
211        end if
212        if (config%do_clear) then
213          flux%sw_dn_clear(jcol,:) = 0.0_jprb
214          flux%sw_up_clear(jcol,:) = 0.0_jprb
215          if (allocated(flux%sw_dn_direct_clear)) then
216            flux%sw_dn_direct_clear(jcol,:) = 0.0_jprb
217          end if
218        end if
219
220        if (config%do_save_spectral_flux) then
221          flux%sw_dn_band(:,jcol,:) = 0.0_jprb
222          flux%sw_up_band(:,jcol,:) = 0.0_jprb
223          if (allocated(flux%sw_dn_direct_band)) then
224            flux%sw_dn_direct_band(:,jcol,:) = 0.0_jprb
225          end if
226          if (config%do_clear) then
227            flux%sw_dn_clear_band(:,jcol,:) = 0.0_jprb
228            flux%sw_up_clear_band(:,jcol,:) = 0.0_jprb
229            if (allocated(flux%sw_dn_direct_clear_band)) then
230              flux%sw_dn_direct_clear_band(:,jcol,:) = 0.0_jprb
231            end if
232          end if
233        end if
234
235        flux%sw_dn_diffuse_surf_g(:,jcol) = 0.0_jprb
236        flux%sw_dn_direct_surf_g(:,jcol)  = 0.0_jprb
237        if (config%do_clear) then
238          flux%sw_dn_diffuse_surf_clear_g(:,jcol) = 0.0_jprb
239          flux%sw_dn_direct_surf_clear_g(:,jcol)  = 0.0_jprb
240        end if
241
242        cycle
243      end if ! sun is below the horizon
244
245      ! At this point mu0 >= 1.0e-10
246
247      ! Define which layers contain cloud; assume that
248      ! cloud%crop_cloud_fraction has already been called
249      is_clear_sky_layer = .true.
250      do jlev = 1,nlev
251        if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
252          is_clear_sky_layer(jlev) = .false.
253        end if
254      end do
255
256      ! --------------------------------------------------------
257      ! Section 3: Loop over layers to compute reflectance and transmittance
258      ! --------------------------------------------------------
259      ! In this section the reflectance, transmittance and sources
260      ! are computed for each layer
261      do jlev = 1,nlev ! Start at top-of-atmosphere
262
263        ! Array-wise assignments
264        gamma1 = 0.0_jprb
265        gamma2 = 0.0_jprb
266        gamma3 = 0.0_jprb
267
268        nreg = nregions
269        if (is_clear_sky_layer(jlev)) then
270          nreg = 1
271        end if
272        do jreg = 1,nreg
273          if (jreg == 1) then
274            od_total  =  od(:,jlev,jcol)
275            ssa_total = ssa(:,jlev,jcol)
276            g_total   =   g(:,jlev,jcol)
277          else
278            do jg = 1,ng
279              ! Mapping from g-point to band
280              iband = config%i_band_from_reordered_g_sw(jg)
281              scat_od = od(jg,jlev,jcol)*ssa(jg,jlev,jcol)
282              scat_od_cloud = od_cloud(iband,jlev,jcol) &
283                   &  * ssa_cloud(iband,jlev,jcol) * od_scaling(jreg,jlev,jcol)
284              ! Add scaled cloud optical depth to clear-sky value
285              od_total(jg) = od(jg,jlev,jcol) &
286                   &  + od_cloud(iband,jlev,jcol)*od_scaling(jreg,jlev,jcol)
287              ! Compute single-scattering albedo and asymmetry
288              ! factor of gas-cloud combination
289              ssa_total(jg) = (scat_od+scat_od_cloud) &
290                   &  / od_total(jg)
291              g_total(jg) = (scat_od*g(jg,jlev,jcol) &
292                   &         + scat_od_cloud * g_cloud(iband,jlev,jcol)) &
293                   &      / (scat_od + scat_od_cloud)
294            end do
295          end if
296
297          if (config%do_sw_delta_scaling_with_gases) then
298            ! Apply delta-Eddington scaling to the aerosol-gas(-cloud)
299            ! mixture
300            call delta_eddington(od_total, ssa_total, g_total)
301          end if
302          call calc_two_stream_gammas_sw(ng, &
303               &  mu0, ssa_total, g_total, &
304               &  gamma1, gamma2, gamma3)
305          call calc_reflectance_transmittance_sw(ng, &
306               &  mu0, od_total, ssa_total, &
307               &  gamma1, gamma2, gamma3, &
308               &  reflectance(:,jreg,jlev), transmittance(:,jreg,jlev), &
309               &  ref_dir(:,jreg,jlev), trans_dir_diff(:,jreg,jlev), &
310               &  trans_dir_dir(:,jreg,jlev) )
311        end do
312      end do
313
314      ! --------------------------------------------------------
315      ! Section 4: Compute total albedos
316      ! --------------------------------------------------------
317
318      total_albedo(:,:,:) = 0.0_jprb
319      total_albedo_direct(:,:,:) = 0.0_jprb
320
321      ! Copy surface albedo in clear-sky region
322      do jg = 1,ng
323        total_albedo(jg,1,nlev+1) = albedo_diffuse(jg,jcol)
324      end do
325
326      ! If direct albedo is available, use it; otherwise copy from
327      ! diffuse albedo
328      do jg = 1,ng
329        total_albedo_direct(jg,1,nlev+1) &
330             &  = mu0 * albedo_direct(jg,jcol)
331      end do
332
333      ! If there is cloud in the lowest layer then we need the albedos
334      ! underneath
335      if (.not. is_clear_sky_layer(nlev)) then
336        do jreg = 2,nregions
337          total_albedo(:,jreg,nlev+1)        = total_albedo(:,1,nlev+1)
338          total_albedo_direct(:,jreg,nlev+1) = total_albedo_direct(:,1,nlev+1)
339        end do
340      end if
341
342      if (config%do_clear) then
343        total_albedo_clear(:,nlev+1)        = total_albedo(:,1,nlev+1)
344        total_albedo_clear_direct(:,nlev+1) = total_albedo_direct(:,1,nlev+1)
345      end if
346
347      ! Work up from the surface computing the total albedo of the
348      ! atmosphere below that point using the adding method
349      do jlev = nlev,1,-1
350
351        total_albedo_below        = 0.0_jprb
352        total_albedo_below_direct = 0.0_jprb
353
354        if (config%do_clear) then
355          ! For clear-skies there is no need to consider "above" and
356          ! "below" quantities since with no cloud overlap to worry
357          ! about, these are the same
358          inv_denom(:,1) = 1.0_jprb &
359               &  / (1.0_jprb - total_albedo_clear(:,jlev+1)*reflectance(:,1,jlev))
360          total_albedo_clear(:,jlev) = reflectance(:,1,jlev) &
361               &  + transmittance(:,1,jlev) * transmittance(:,1,jlev) &
362               &  * total_albedo_clear(:,jlev+1) * inv_denom(:,1)
363          total_albedo_clear_direct(:,jlev) = ref_dir(:,1,jlev) &
364               &  + (trans_dir_dir(:,1,jlev)*total_albedo_clear_direct(:,jlev+1) &
365               &     +trans_dir_diff(:,1,jlev)*total_albedo_clear(:,jlev+1)) &
366               &  * transmittance(:,1,jlev) * inv_denom(:,1)
367        end if
368
369        if (is_clear_sky_layer(jlev)) then
370          inv_denom(:,1) = 1.0_jprb &
371               &  / (1.0_jprb - total_albedo(:,1,jlev+1)*reflectance(:,1,jlev))
372          total_albedo_below(:,1) = reflectance(:,1,jlev) &
373               &  + transmittance(:,1,jlev)  * transmittance(:,1,jlev) &
374               &  * total_albedo(:,1,jlev+1) * inv_denom(:,1)
375          total_albedo_below_direct(:,1) = ref_dir(:,1,jlev) &
376               &  + (trans_dir_dir(:,1,jlev)*total_albedo_direct(:,1,jlev+1) &
377               &     +trans_dir_diff(:,1,jlev)*total_albedo(:,1,jlev+1)) &
378               &  * transmittance(:,1,jlev) * inv_denom(:,1)
379        else
380          inv_denom = 1.0_jprb / (1.0_jprb - total_albedo(:,:,jlev+1)*reflectance(:,:,jlev))
381          total_albedo_below = reflectance(:,:,jlev) &
382               &  + transmittance(:,:,jlev)  * transmittance(:,:,jlev) &
383               &  * total_albedo(:,:,jlev+1) * inv_denom
384          total_albedo_below_direct = ref_dir(:,:,jlev) &
385               &  + (trans_dir_dir(:,:,jlev)*total_albedo_direct(:,:,jlev+1) &
386               &     +trans_dir_diff(:,:,jlev)*total_albedo(:,:,jlev+1)) &
387               &  * transmittance(:,:,jlev) * inv_denom
388        end if
389
390        ! Account for cloud overlap when converting albedo below a
391        ! layer interface to the equivalent values just above
392        if (is_clear_sky_layer(jlev) .and. is_clear_sky_layer(jlev-1)) then
393          total_albedo(:,:,jlev)        = total_albedo_below(:,:)
394          total_albedo_direct(:,:,jlev) = total_albedo_below_direct(:,:)
395        else
396          ! Use overlap matrix and exclude "anomalous" horizontal
397          ! transport described by Shonk & Hogan (2008).  Therefore,
398          ! the operation we perform is essentially diag(total_albedo)
399          ! = matmul(transpose(v_matrix)), diag(total_albedo_below)).
400          do jreg = 1,nregions
401            do jreg2 = 1,nregions
402              total_albedo(:,jreg,jlev) &
403                   &  = total_albedo(:,jreg,jlev) &
404                   &  + total_albedo_below(:,jreg2) &
405                   &  * v_matrix(jreg2,jreg,jlev,jcol)
406              total_albedo_direct(:,jreg,jlev) &
407                   &  = total_albedo_direct(:,jreg,jlev) &
408                   &  + total_albedo_below_direct(:,jreg2) &
409                   &  * v_matrix(jreg2,jreg,jlev,jcol)
410            end do
411          end do
412
413        end if
414
415      end do ! Reverse loop over levels
416
417      ! --------------------------------------------------------
418      ! Section 5: Compute fluxes
419      ! --------------------------------------------------------
420
421      ! Top-of-atmosphere fluxes into the regions of the top-most
422      ! layer, zero since we assume no diffuse downwelling
423      flux_dn = 0.0_jprb
424      ! Direct downwelling flux (into a plane perpendicular to the
425      ! sun) entering the top of each region in the top-most layer
426      do jreg = 1,nregions
427        direct_dn(:,jreg) = incoming_sw(:,jcol)*region_fracs(jreg,1,jcol)
428      end do
429      flux_up = direct_dn*total_albedo_direct(:,:,1)
430
431      if (config%do_clear) then
432        flux_dn_clear = 0.0_jprb
433        direct_dn_clear(:) = incoming_sw(:,jcol)
434        flux_up_clear = direct_dn_clear*total_albedo_clear_direct(:,1)
435      end if
436
437      ! Store the TOA broadband fluxes
438      flux%sw_up(jcol,1) = sum(sum(flux_up,1))
439      flux%sw_dn(jcol,1) = mu0 * sum(sum(direct_dn,1))
440      if (allocated(flux%sw_dn_direct)) then
441        flux%sw_dn_direct(jcol,1) = flux%sw_dn(jcol,1)
442      end if
443      if (config%do_clear) then
444        flux%sw_up_clear(jcol,1) = sum(flux_up_clear)
445        flux%sw_dn_clear(jcol,1) = mu0 * sum(direct_dn_clear)
446        if (allocated(flux%sw_dn_direct_clear)) then
447          flux%sw_dn_direct_clear(jcol,1) = flux%sw_dn_clear(jcol,1)
448        end if
449      end if
450
451      ! Save the spectral fluxes if required; some redundancy here as
452      ! the TOA downwelling flux is the same in clear and cloudy skies
453      if (config%do_save_spectral_flux) then
454        call indexed_sum(sum(flux_up,2), &
455             &           config%i_spec_from_reordered_g_sw, &
456             &           flux%sw_up_band(:,jcol,1))
457        call indexed_sum(sum(direct_dn,2), &
458             &           config%i_spec_from_reordered_g_sw, &
459             &           flux%sw_dn_band(:,jcol,1))
460        flux%sw_dn_band(:,jcol,1) = &
461             &  mu0 * flux%sw_dn_band(:,jcol,1)
462        if (allocated(flux%sw_dn_direct_band)) then
463          flux%sw_dn_direct_band(:,jcol,1) = flux%sw_dn_band(:,jcol,1)
464        end if
465        call add_indexed_sum(sum(flux_dn,2), &
466             &           config%i_spec_from_reordered_g_sw, &
467             &           flux%sw_dn_band(:,jcol,1))
468        if (config%do_clear) then
469          call indexed_sum(flux_up_clear, &
470               &           config%i_spec_from_reordered_g_sw, &
471               &           flux%sw_up_clear_band(:,jcol,1))
472          call indexed_sum(direct_dn_clear, &
473               &           config%i_spec_from_reordered_g_sw, &
474               &           flux%sw_dn_clear_band(:,jcol,1))
475          flux%sw_dn_clear_band(:,jcol,1) &
476               &  = mu0 * flux%sw_dn_clear_band(:,jcol,1)
477          if (allocated(flux%sw_dn_direct_clear_band)) then
478            flux%sw_dn_direct_clear_band(:,jcol,1) &
479                 &  = flux%sw_dn_clear_band(:,jcol,1)
480          end if
481          call add_indexed_sum(flux_dn_clear, &
482               &           config%i_spec_from_reordered_g_sw, &
483               &           flux%sw_dn_clear_band(:,jcol,1))
484        end if
485      end if
486
487      ! Final loop back down through the atmosphere to compute fluxes
488      do jlev = 1,nlev
489        if (config%do_clear) then
490          flux_dn_clear = (transmittance(:,1,jlev)*flux_dn_clear + direct_dn_clear &
491               &  * (trans_dir_dir(:,1,jlev)*total_albedo_clear_direct(:,jlev+1)*reflectance(:,1,jlev) &
492               &     + trans_dir_diff(:,1,jlev) )) &
493               &  / (1.0_jprb - reflectance(:,1,jlev)*total_albedo_clear(:,jlev+1))
494          direct_dn_clear = trans_dir_dir(:,1,jlev)*direct_dn_clear
495          flux_up_clear = direct_dn_clear*total_albedo_clear_direct(:,jlev+1) &
496               &        +   flux_dn_clear*total_albedo_clear(:,jlev+1)
497        end if
498
499        if (is_clear_sky_layer(jlev)) then
500          flux_dn(:,1) = (transmittance(:,1,jlev)*flux_dn(:,1) + direct_dn(:,1) &
501               &  * (trans_dir_dir(:,1,jlev)*total_albedo_direct(:,1,jlev+1)*reflectance(:,1,jlev) &
502               &     + trans_dir_diff(:,1,jlev) )) &
503               &  / (1.0_jprb - reflectance(:,1,jlev)*total_albedo(:,1,jlev+1))
504          direct_dn(:,1) = trans_dir_dir(:,1,jlev)*direct_dn(:,1)
505          flux_up(:,1) = direct_dn(:,1)*total_albedo_direct(:,1,jlev+1) &
506               &  +        flux_dn(:,1)*total_albedo(:,1,jlev+1)
507          flux_dn(:,2:)  = 0.0_jprb
508          flux_up(:,2:)  = 0.0_jprb
509          direct_dn(:,2:)= 0.0_jprb
510        else
511          flux_dn = (transmittance(:,:,jlev)*flux_dn + direct_dn &
512               &  * (trans_dir_dir(:,:,jlev)*total_albedo_direct(:,:,jlev+1)*reflectance(:,:,jlev) &
513               &     + trans_dir_diff(:,:,jlev) )) &
514               &  / (1.0_jprb - reflectance(:,:,jlev)*total_albedo(:,:,jlev+1))
515          direct_dn = trans_dir_dir(:,:,jlev)*direct_dn
516          flux_up = direct_dn*total_albedo_direct(:,:,jlev+1) &
517               &  +   flux_dn*total_albedo(:,:,jlev+1)
518        end if
519
520        if (.not. (is_clear_sky_layer(jlev) &
521             &    .and. is_clear_sky_layer(jlev+1))) then
522          ! Account for overlap rules in translating fluxes just above
523          ! a layer interface to the values just below
524          flux_dn_below = singlemat_x_vec(ng,ng,nregions, &
525               &  v_matrix(:,:,jlev+1,jcol), flux_dn)
526          direct_dn_below = singlemat_x_vec(ng,ng,nregions, &
527               &  v_matrix(:,:,jlev+1,jcol), direct_dn)
528          flux_dn = flux_dn_below
529          direct_dn = direct_dn_below
530        end if ! Otherwise the fluxes in each region are the same so
531               ! nothing to do
532
533        ! Store the broadband fluxes
534        flux%sw_up(jcol,jlev+1) = sum(sum(flux_up,1))
535        if (allocated(flux%sw_dn_direct)) then
536          flux%sw_dn_direct(jcol,jlev+1) = mu0 * sum(sum(direct_dn,1))
537          flux%sw_dn(jcol,jlev+1) &
538               &  = flux%sw_dn_direct(jcol,jlev+1) + sum(sum(flux_dn,1))
539        else
540          flux%sw_dn(jcol,jlev+1) = mu0 * sum(sum(direct_dn,1)) + sum(sum(flux_dn,1))   
541        end if
542        if (config%do_clear) then
543          flux%sw_up_clear(jcol,jlev+1) = sum(flux_up_clear)
544          if (allocated(flux%sw_dn_direct_clear)) then
545            flux%sw_dn_direct_clear(jcol,jlev+1) = mu0 * sum(direct_dn_clear)
546            flux%sw_dn_clear(jcol,jlev+1) &
547                 &  = flux%sw_dn_direct_clear(jcol,jlev+1) + sum(flux_dn_clear)
548          else
549            flux%sw_dn_clear(jcol,jlev+1) = mu0 * sum(direct_dn_clear) &
550                 &  + sum(flux_dn_clear)
551          end if
552        end if
553
554        ! Save the spectral fluxes if required
555        if (config%do_save_spectral_flux) then
556          call indexed_sum(sum(flux_up,2), &
557               &           config%i_spec_from_reordered_g_sw, &
558               &           flux%sw_up_band(:,jcol,jlev+1))
559          call indexed_sum(sum(direct_dn,2), &
560               &           config%i_spec_from_reordered_g_sw, &
561               &           flux%sw_dn_band(:,jcol,jlev+1))
562          flux%sw_dn_band(:,jcol,jlev+1) = &
563               &  mu0 * flux%sw_dn_band(:,jcol,jlev+1)
564          if (allocated(flux%sw_dn_direct_band)) then
565            flux%sw_dn_direct_band(:,jcol,jlev+1) &
566                 &  = flux%sw_dn_band(:,jcol,jlev+1)
567          end if
568          call add_indexed_sum(sum(flux_dn,2), &
569               &           config%i_spec_from_reordered_g_sw, &
570               &           flux%sw_dn_band(:,jcol,jlev+1))
571          if (config%do_clear) then
572            call indexed_sum(flux_up_clear, &
573                 &           config%i_spec_from_reordered_g_sw, &
574                 &           flux%sw_up_clear_band(:,jcol,jlev+1))
575            call indexed_sum(direct_dn_clear, &
576                 &           config%i_spec_from_reordered_g_sw, &
577                 &           flux%sw_dn_clear_band(:,jcol,jlev+1))
578            flux%sw_dn_clear_band(:,jcol,jlev+1) = &
579                 &  mu0 * flux%sw_dn_clear_band(:,jcol,jlev+1)
580            if (allocated(flux%sw_dn_direct_clear_band)) then
581              flux%sw_dn_direct_clear_band(:,jcol,jlev+1) &
582                   &  = flux%sw_dn_clear_band(:,jcol,jlev+1)
583            end if
584            call add_indexed_sum(flux_dn_clear, &
585                 &           config%i_spec_from_reordered_g_sw, &
586                 &           flux%sw_dn_clear_band(:,jcol,jlev+1))
587          end if
588        end if
589
590      end do ! Final loop over levels
591     
592      ! Store surface spectral fluxes, if required (after the end of
593      ! the final loop over levels, the current values of these arrays
594      ! will be the surface values)
595      flux%sw_dn_diffuse_surf_g(:,jcol) = sum(flux_dn,2)
596      flux%sw_dn_direct_surf_g(:,jcol)  = mu0 * sum(direct_dn,2)
597      if (config%do_clear) then
598        flux%sw_dn_diffuse_surf_clear_g(:,jcol) = flux_dn_clear
599        flux%sw_dn_direct_surf_clear_g(:,jcol)  = mu0 * direct_dn_clear
600      end if
601
602    end do ! Loop over columns
603
604    if (lhook) call dr_hook('radiation_tripleclouds_sw:solver_tripleclouds_sw',1,hook_handle)
605
606  end subroutine solver_tripleclouds_sw
607
608end module radiation_tripleclouds_sw
Note: See TracBrowser for help on using the repository browser.