source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_mcica_sw.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: 14.4 KB
Line 
1! radiation_mcica_sw.F90 - Monte-Carlo Independent Column Approximation shortwave solver
2!
3! (C) Copyright 2015- 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
20module radiation_mcica_sw
21
22  public
23
24contains
25
26  ! Provides elemental function "delta_eddington"
27#include "radiation_delta_eddington.h"
28
29  !---------------------------------------------------------------------
30  ! Shortwave Monte Carlo Independent Column Approximation
31  ! (McICA). This implementation performs a clear-sky and a cloudy-sky
32  ! calculation, and then weights the two to get the all-sky fluxes
33  ! according to the total cloud cover. This method reduces noise for
34  ! low cloud cover situations, and exploits the clear-sky
35  ! calculations that are usually performed for diagnostic purposes
36  ! simultaneously. The cloud generator has been carefully written
37  ! such that the stochastic cloud field satisfies the prescribed
38  ! overlap parameter accounting for this weighting.
39  subroutine solver_mcica_sw(nlev,istartcol,iendcol, &
40       &  config, single_level, cloud, &
41       &  od, ssa, g, od_cloud, ssa_cloud, g_cloud, &
42       &  albedo_direct, albedo_diffuse, incoming_sw, &
43       &  flux)
44
45    use parkind1, only           : jprb
46    use yomhook,  only           : lhook, dr_hook, jphook
47
48    use radiation_io,   only           : nulerr, radiation_abort
49    use radiation_config, only         : config_type
50    use radiation_single_level, only   : single_level_type
51    use radiation_cloud, only          : cloud_type
52    use radiation_flux, only           : flux_type
53    use radiation_two_stream, only     : calc_ref_trans_sw
54    use radiation_adding_ica_sw, only  : adding_ica_sw
55    use radiation_cloud_generator, only: cloud_generator
56
57    implicit none
58
59    ! Inputs
60    integer, intent(in) :: nlev               ! number of model levels
61    integer, intent(in) :: istartcol, iendcol ! range of columns to process
62    type(config_type),        intent(in) :: config
63    type(single_level_type),  intent(in) :: single_level
64    type(cloud_type),         intent(in) :: cloud
65
66    ! Gas and aerosol optical depth, single-scattering albedo and
67    ! asymmetry factor at each shortwave g-point
68    real(jprb), intent(in), dimension(config%n_g_sw, nlev, istartcol:iendcol) :: &
69         &  od, ssa, g
70
71    ! Cloud and precipitation optical depth, single-scattering albedo and
72    ! asymmetry factor in each shortwave band
73    real(jprb), intent(in), dimension(config%n_bands_sw,nlev,istartcol:iendcol)   :: &
74         &  od_cloud, ssa_cloud, g_cloud
75
76    ! Direct and diffuse surface albedos, and the incoming shortwave
77    ! flux into a plane perpendicular to the incoming radiation at
78    ! top-of-atmosphere in each of the shortwave g points
79    real(jprb), intent(in), dimension(config%n_g_sw,istartcol:iendcol) :: &
80         &  albedo_direct, albedo_diffuse, incoming_sw
81
82    ! Output
83    type(flux_type), intent(inout):: flux
84
85    ! Local variables
86
87    ! Cosine of solar zenith angle
88    real(jprb)                                 :: cos_sza
89
90    ! Diffuse reflectance and transmittance for each layer in clear
91    ! and all skies
92    real(jprb), dimension(config%n_g_sw, nlev) :: ref_clear, trans_clear, reflectance, transmittance
93
94    ! Fraction of direct beam scattered by a layer into the upwelling
95    ! or downwelling diffuse streams, in clear and all skies
96    real(jprb), dimension(config%n_g_sw, nlev) :: ref_dir_clear, trans_dir_diff_clear, ref_dir, trans_dir_diff
97
98    ! Transmittance for the direct beam in clear and all skies
99    real(jprb), dimension(config%n_g_sw, nlev) :: trans_dir_dir_clear, trans_dir_dir
100
101    ! Fluxes per g point
102    real(jprb), dimension(config%n_g_sw, nlev+1) :: flux_up, flux_dn_diffuse, flux_dn_direct
103
104    ! Combined gas+aerosol+cloud optical depth, single scattering
105    ! albedo and asymmetry factor
106    real(jprb), dimension(config%n_g_sw) :: od_total, ssa_total, g_total
107
108    ! Combined scattering optical depth
109    real(jprb) :: scat_od
110
111    ! Optical depth scaling from the cloud generator, zero indicating
112    ! clear skies
113    real(jprb), dimension(config%n_g_sw,nlev) :: od_scaling
114
115    ! Modified optical depth after McICA scaling to represent cloud
116    ! inhomogeneity
117    real(jprb), dimension(config%n_g_sw) :: od_cloud_new
118
119    ! Total cloud cover output from the cloud generator
120    real(jprb) :: total_cloud_cover
121
122    ! Number of g points
123    integer :: ng
124
125    ! Loop indices for level, column and g point
126    integer :: jlev, jcol, jg
127
128    real(jphook) :: hook_handle
129
130    if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',0,hook_handle)
131
132    if (.not. config%do_clear) then
133      write(nulerr,'(a)') '*** Error: shortwave McICA requires clear-sky calculation to be performed'
134      call radiation_abort()     
135    end if
136
137    ng = config%n_g_sw
138
139    ! Loop through columns
140    do jcol = istartcol,iendcol
141      ! Only perform calculation if sun above the horizon
142      if (single_level%cos_sza(jcol) > 0.0_jprb) then
143        cos_sza = single_level%cos_sza(jcol)
144
145        ! Clear-sky calculation - first compute clear-sky reflectance,
146        ! transmittance etc at each model level
147        if (.not. config%do_sw_delta_scaling_with_gases) then
148          ! Delta-Eddington scaling has already been performed to the
149          ! aerosol part of od, ssa and g
150          call calc_ref_trans_sw(ng*nlev, &
151               &  cos_sza, od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), &
152               &  ref_clear, trans_clear, &
153               &  ref_dir_clear, trans_dir_diff_clear, &
154               &  trans_dir_dir_clear)
155        else
156          ! Apply delta-Eddington scaling to the aerosol-gas mixture
157          do jlev = 1,nlev
158            od_total  =  od(:,jlev,jcol)
159            ssa_total = ssa(:,jlev,jcol)
160            g_total   =   g(:,jlev,jcol)
161            call delta_eddington(od_total, ssa_total, g_total)
162            call calc_ref_trans_sw(ng, &
163                 &  cos_sza, od_total, ssa_total, g_total, &
164                 &  ref_clear(:,jlev), trans_clear(:,jlev), &
165                 &  ref_dir_clear(:,jlev), trans_dir_diff_clear(:,jlev), &
166                 &  trans_dir_dir_clear(:,jlev) )
167          end do
168        end if
169
170        ! Use adding method to compute fluxes
171        call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), &
172             &  albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), &
173             &  ref_clear, trans_clear, ref_dir_clear, trans_dir_diff_clear, &
174             &  trans_dir_dir_clear, flux_up, flux_dn_diffuse, flux_dn_direct)
175       
176        ! Sum over g-points to compute and save clear-sky broadband
177        ! fluxes
178        flux%sw_up_clear(jcol,:) = sum(flux_up,1)
179        if (allocated(flux%sw_dn_direct_clear)) then
180          flux%sw_dn_direct_clear(jcol,:) &
181               &  = sum(flux_dn_direct,1)
182          flux%sw_dn_clear(jcol,:) = sum(flux_dn_diffuse,1) &
183               &  + flux%sw_dn_direct_clear(jcol,:)
184        else
185          flux%sw_dn_clear(jcol,:) = sum(flux_dn_diffuse,1) &
186               &  + sum(flux_dn_direct,1)
187        end if
188        ! Store spectral downwelling fluxes at surface
189        flux%sw_dn_diffuse_surf_clear_g(:,jcol) = flux_dn_diffuse(:,nlev+1)
190        flux%sw_dn_direct_surf_clear_g(:,jcol)  = flux_dn_direct(:,nlev+1)
191
192        ! Do cloudy-sky calculation
193        call cloud_generator(ng, nlev, config%i_overlap_scheme, &
194             &  single_level%iseed(jcol), &
195             &  config%cloud_fraction_threshold, &
196             &  cloud%fraction(jcol,:), cloud%overlap_param(jcol,:), &
197             &  config%cloud_inhom_decorr_scaling, cloud%fractional_std(jcol,:), &
198             &  config%pdf_sampler, od_scaling, total_cloud_cover, &
199             &  use_beta_overlap=config%use_beta_overlap, &
200             &  use_vectorizable_generator=config%use_vectorizable_generator)
201
202        ! Store total cloud cover
203        flux%cloud_cover_sw(jcol) = total_cloud_cover
204       
205        if (total_cloud_cover >= config%cloud_fraction_threshold) then
206          ! Total-sky calculation
207          do jlev = 1,nlev
208            ! Compute combined gas+aerosol+cloud optical properties
209            if (cloud%fraction(jcol,jlev) >= config%cloud_fraction_threshold) then
210              do jg = 1,ng
211                od_cloud_new(jg) = od_scaling(jg,jlev) &
212                   &  * od_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol)
213                od_total(jg)  = od(jg,jlev,jcol) + od_cloud_new(jg)
214                ssa_total(jg) = 0.0_jprb
215                g_total(jg)   = 0.0_jprb
216
217                ! In single precision we need to protect against the
218                ! case that od_total > 0.0 and ssa_total > 0.0 but
219                ! od_total*ssa_total == 0 due to underflow
220                if (od_total(jg) > 0.0_jprb) then
221                  scat_od = ssa(jg,jlev,jcol)*od(jg,jlev,jcol) &
222                       &     + ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
223                       &     *  od_cloud_new(jg)
224                  ssa_total(jg) = scat_od / od_total(jg)
225                  if (scat_od > 0.0_jprb) then
226                    g_total(jg) = (g(jg,jlev,jcol)*ssa(jg,jlev,jcol)*od(jg,jlev,jcol) &
227                         &     +   g_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
228                         &     * ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
229                         &     *  od_cloud_new(jg)) &
230                         &     / scat_od
231                  end if
232                end if
233              end do
234
235              ! Apply delta-Eddington scaling to the cloud-aerosol-gas
236              ! mixture
237              if (config%do_sw_delta_scaling_with_gases) then
238                call delta_eddington(od_total, ssa_total, g_total)
239              end if
240
241              ! Compute cloudy-sky reflectance, transmittance etc at
242              ! each model level
243              call calc_ref_trans_sw(ng, &
244                   &  cos_sza, od_total, ssa_total, g_total, &
245                   &  reflectance(:,jlev), transmittance(:,jlev), &
246                   &  ref_dir(:,jlev), trans_dir_diff(:,jlev), &
247                   &  trans_dir_dir(:,jlev))
248             
249            else
250              ! Clear-sky layer: copy over clear-sky values
251              reflectance(:,jlev) = ref_clear(:,jlev)
252              transmittance(:,jlev) = trans_clear(:,jlev)
253              ref_dir(:,jlev) = ref_dir_clear(:,jlev)
254              trans_dir_diff(:,jlev) = trans_dir_diff_clear(:,jlev)
255              trans_dir_dir(:,jlev) = trans_dir_dir_clear(:,jlev)
256            end if
257          end do
258           
259          ! Use adding method to compute fluxes for an overcast sky
260          call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), &
261               &  albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), &
262               &  reflectance, transmittance, ref_dir, trans_dir_diff, &
263               &  trans_dir_dir, flux_up, flux_dn_diffuse, flux_dn_direct)
264         
265          ! Store overcast broadband fluxes
266          flux%sw_up(jcol,:) = sum(flux_up,1)
267          if (allocated(flux%sw_dn_direct)) then
268            flux%sw_dn_direct(jcol,:) = sum(flux_dn_direct,1)
269            flux%sw_dn(jcol,:) = sum(flux_dn_diffuse,1) &
270                 &  + flux%sw_dn_direct(jcol,:)
271          else
272            flux%sw_dn(jcol,:) = sum(flux_dn_diffuse,1) &
273                 &  + sum(flux_dn_direct,1)
274          end if
275
276          ! Cloudy flux profiles currently assume completely overcast
277          ! skies; perform weighted average with clear-sky profile
278          flux%sw_up(jcol,:) =  total_cloud_cover *flux%sw_up(jcol,:) &
279               &  + (1.0_jprb - total_cloud_cover)*flux%sw_up_clear(jcol,:)
280          flux%sw_dn(jcol,:) =  total_cloud_cover *flux%sw_dn(jcol,:) &
281               &  + (1.0_jprb - total_cloud_cover)*flux%sw_dn_clear(jcol,:)
282          if (allocated(flux%sw_dn_direct)) then
283            flux%sw_dn_direct(jcol,:) = total_cloud_cover *flux%sw_dn_direct(jcol,:) &
284                 &  + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_clear(jcol,:)
285          end if
286          ! Likewise for surface spectral fluxes
287          flux%sw_dn_diffuse_surf_g(:,jcol) = flux_dn_diffuse(:,nlev+1)
288          flux%sw_dn_direct_surf_g(:,jcol)  = flux_dn_direct(:,nlev+1)
289          flux%sw_dn_diffuse_surf_g(:,jcol) = total_cloud_cover *flux%sw_dn_diffuse_surf_g(:,jcol) &
290               &     + (1.0_jprb - total_cloud_cover)*flux%sw_dn_diffuse_surf_clear_g(:,jcol)
291          flux%sw_dn_direct_surf_g(:,jcol) = total_cloud_cover *flux%sw_dn_direct_surf_g(:,jcol) &
292               &     + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_surf_clear_g(:,jcol)
293         
294        else
295          ! No cloud in profile and clear-sky fluxes already
296          ! calculated: copy them over
297          flux%sw_up(jcol,:) = flux%sw_up_clear(jcol,:)
298          flux%sw_dn(jcol,:) = flux%sw_dn_clear(jcol,:)
299          if (allocated(flux%sw_dn_direct)) then
300            flux%sw_dn_direct(jcol,:) = flux%sw_dn_direct_clear(jcol,:)
301          end if
302          flux%sw_dn_diffuse_surf_g(:,jcol) = flux%sw_dn_diffuse_surf_clear_g(:,jcol)
303          flux%sw_dn_direct_surf_g(:,jcol)  = flux%sw_dn_direct_surf_clear_g(:,jcol)
304
305        end if ! Cloud is present in profile
306
307      else
308        ! Set fluxes to zero if sun is below the horizon
309        flux%sw_up(jcol,:) = 0.0_jprb
310        flux%sw_dn(jcol,:) = 0.0_jprb
311        if (allocated(flux%sw_dn_direct)) then
312          flux%sw_dn_direct(jcol,:) = 0.0_jprb
313        end if
314        flux%sw_up_clear(jcol,:) = 0.0_jprb
315        flux%sw_dn_clear(jcol,:) = 0.0_jprb
316        if (allocated(flux%sw_dn_direct_clear)) then
317          flux%sw_dn_direct_clear(jcol,:) = 0.0_jprb
318        end if
319        flux%sw_dn_diffuse_surf_g(:,jcol) = 0.0_jprb
320        flux%sw_dn_direct_surf_g(:,jcol)  = 0.0_jprb
321        flux%sw_dn_diffuse_surf_clear_g(:,jcol) = 0.0_jprb
322        flux%sw_dn_direct_surf_clear_g(:,jcol)  = 0.0_jprb
323      end if ! Sun above horizon
324
325    end do ! Loop over columns
326
327    if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',1,hook_handle)
328   
329  end subroutine solver_mcica_sw
330
331end module radiation_mcica_sw
Note: See TracBrowser for help on using the repository browser.