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

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

Online implementation of the radiative transfer code ECRAD in LMDZ.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
  • Adaptation of compilation scripts (CPP_ECRAD keys)
  • Call of ecrad in radlwsw_m.F90 under the logical key iflag_rrtm = 2
File size: 15.1 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
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_two_stream_gammas_sw, &
54         &                               calc_reflectance_transmittance_sw
55    use radiation_adding_ica_sw, only  : adding_ica_sw
56    use radiation_cloud_generator, only: cloud_generator
57
58    implicit none
59
60    ! Inputs
61    integer, intent(in) :: nlev               ! number of model levels
62    integer, intent(in) :: istartcol, iendcol ! range of columns to process
63    type(config_type),        intent(in) :: config
64    type(single_level_type),  intent(in) :: single_level
65    type(cloud_type),         intent(in) :: cloud
66
67    ! Gas and aerosol optical depth, single-scattering albedo and
68    ! asymmetry factor at each shortwave g-point
69    real(jprb), intent(in), dimension(config%n_g_sw, nlev, istartcol:iendcol) :: &
70         &  od, ssa, g
71
72    ! Cloud and precipitation optical depth, single-scattering albedo and
73    ! asymmetry factor in each shortwave band
74    real(jprb), intent(in), dimension(config%n_bands_sw,nlev,istartcol:iendcol)   :: &
75         &  od_cloud, ssa_cloud, g_cloud
76
77    ! Direct and diffuse surface albedos, and the incoming shortwave
78    ! flux into a plane perpendicular to the incoming radiation at
79    ! top-of-atmosphere in each of the shortwave g points
80    real(jprb), intent(in), dimension(config%n_g_sw,istartcol:iendcol) :: &
81         &  albedo_direct, albedo_diffuse, incoming_sw
82
83    ! Output
84    type(flux_type), intent(inout):: flux
85
86    ! Local variables
87
88    ! Cosine of solar zenith angle
89    real(jprb)                                 :: cos_sza
90
91    ! Diffuse reflectance and transmittance for each layer in clear
92    ! and all skies
93    real(jprb), dimension(config%n_g_sw, nlev) :: ref_clear, trans_clear, reflectance, transmittance
94
95    ! Fraction of direct beam scattered by a layer into the upwelling
96    ! or downwelling diffuse streams, in clear and all skies
97    real(jprb), dimension(config%n_g_sw, nlev) :: ref_dir_clear, trans_dir_diff_clear, ref_dir, trans_dir_diff
98
99    ! Transmittance for the direct beam in clear and all skies
100    real(jprb), dimension(config%n_g_sw, nlev) :: trans_dir_dir_clear, trans_dir_dir
101
102    ! Fluxes per g point
103    real(jprb), dimension(config%n_g_sw, nlev+1) :: flux_up, flux_dn_diffuse, flux_dn_direct
104
105    ! Combined gas+aerosol+cloud optical depth, single scattering
106    ! albedo and asymmetry factor
107    real(jprb), dimension(config%n_g_sw) :: od_total, ssa_total, g_total
108
109    ! Combined scattering optical depth
110    real(jprb) :: scat_od
111
112    ! Two-stream coefficients
113    real(jprb), dimension(config%n_g_sw) :: gamma1, gamma2, gamma3
114
115    ! Optical depth scaling from the cloud generator, zero indicating
116    ! clear skies
117    real(jprb), dimension(config%n_g_sw,nlev) :: od_scaling
118
119    ! Modified optical depth after McICA scaling to represent cloud
120    ! inhomogeneity
121    real(jprb), dimension(config%n_g_sw) :: od_cloud_new
122
123    ! Total cloud cover output from the cloud generator
124    real(jprb) :: total_cloud_cover
125
126    ! Number of g points
127    integer :: ng
128
129    ! Loop indices for level, column and g point
130    integer :: jlev, jcol, jg
131
132    real(jprb) :: hook_handle
133
134    if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',0,hook_handle)
135
136    if (.not. config%do_clear) then
137      write(nulerr,'(a)') '*** Error: shortwave McICA requires clear-sky calculation to be performed'
138      call radiation_abort()     
139    end if
140
141    ng = config%n_g_sw
142
143    ! Loop through columns
144    do jcol = istartcol,iendcol
145      ! Only perform calculation if sun above the horizon
146      if (single_level%cos_sza(jcol) > 0.0_jprb) then
147        cos_sza = single_level%cos_sza(jcol)
148
149        ! Clear-sky calculation - first compute clear-sky reflectance,
150        ! transmittance etc at each model level
151        if (.not. config%do_sw_delta_scaling_with_gases) then
152          ! Delta-Eddington scaling has already been performed to the
153          ! aerosol part of od, ssa and g
154          do jlev = 1,nlev
155            call calc_two_stream_gammas_sw(ng, &
156                 &  cos_sza, ssa(:,jlev,jcol), g(:,jlev,jcol), &
157                 &  gamma1, gamma2, gamma3)
158            call calc_reflectance_transmittance_sw(ng, &
159                 &  cos_sza, od(:,jlev,jcol), ssa(:,jlev,jcol), &
160                 &  gamma1, gamma2, gamma3, &
161                 &  ref_clear(:,jlev), trans_clear(:,jlev), &
162                 &  ref_dir_clear(:,jlev), trans_dir_diff_clear(:,jlev), &
163                 &  trans_dir_dir_clear(:,jlev) )
164          end do
165        else
166          ! Apply delta-Eddington scaling to the aerosol-gas mixture
167          do jlev = 1,nlev
168            od_total  =  od(:,jlev,jcol)
169            ssa_total = ssa(:,jlev,jcol)
170            g_total   =   g(:,jlev,jcol)
171            call delta_eddington(od_total, ssa_total, g_total)
172            call calc_two_stream_gammas_sw(ng, &
173                 &  cos_sza, ssa_total, g_total, &
174                 &  gamma1, gamma2, gamma3)
175            call calc_reflectance_transmittance_sw(ng, &
176                 &  cos_sza, od_total, ssa_total, &
177                 &  gamma1, gamma2, gamma3, &
178                 &  ref_clear(:,jlev), trans_clear(:,jlev), &
179                 &  ref_dir_clear(:,jlev), trans_dir_diff_clear(:,jlev), &
180                 &  trans_dir_dir_clear(:,jlev) )
181          end do
182        end if
183
184        ! Use adding method to compute fluxes
185        call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), &
186             &  albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), &
187             &  ref_clear, trans_clear, ref_dir_clear, trans_dir_diff_clear, &
188             &  trans_dir_dir_clear, flux_up, flux_dn_diffuse, flux_dn_direct)
189       
190        ! Sum over g-points to compute and save clear-sky broadband
191        ! fluxes
192        flux%sw_up_clear(jcol,:) = sum(flux_up,1)
193        if (allocated(flux%sw_dn_direct_clear)) then
194          flux%sw_dn_direct_clear(jcol,:) &
195               &  = sum(flux_dn_direct,1)
196          flux%sw_dn_clear(jcol,:) = sum(flux_dn_diffuse,1) &
197               &  + flux%sw_dn_direct_clear(jcol,:)
198        else
199          flux%sw_dn_clear(jcol,:) = sum(flux_dn_diffuse,1) &
200               &  + sum(flux_dn_direct,1)
201        end if
202        ! Store spectral downwelling fluxes at surface
203        flux%sw_dn_diffuse_surf_clear_g(:,jcol) = flux_dn_diffuse(:,nlev+1)
204        flux%sw_dn_direct_surf_clear_g(:,jcol)  = flux_dn_direct(:,nlev+1)
205
206        ! Do cloudy-sky calculation
207        call cloud_generator(ng, nlev, config%i_overlap_scheme, &
208             &  single_level%iseed(jcol), &
209             &  config%cloud_fraction_threshold, &
210             &  cloud%fraction(jcol,:), cloud%overlap_param(jcol,:), &
211             &  config%cloud_inhom_decorr_scaling, cloud%fractional_std(jcol,:), &
212             &  config%pdf_sampler, od_scaling, total_cloud_cover, &
213             &  is_beta_overlap=config%use_beta_overlap)
214
215        ! Store total cloud cover
216        flux%cloud_cover_sw(jcol) = total_cloud_cover
217       
218        if (total_cloud_cover >= config%cloud_fraction_threshold) then
219          ! Total-sky calculation
220          do jlev = 1,nlev
221            ! Compute combined gas+aerosol+cloud optical properties
222            if (cloud%fraction(jcol,jlev) >= config%cloud_fraction_threshold) then
223              od_cloud_new = od_scaling(:,jlev) &
224                   &  * od_cloud(config%i_band_from_reordered_g_sw,jlev,jcol)
225              od_total  = od(:,jlev,jcol) + od_cloud_new
226              ssa_total = 0.0_jprb
227              g_total   = 0.0_jprb
228              ! In single precision we need to protect against the
229              ! case that od_total > 0.0 and ssa_total > 0.0 but
230              ! od_total*ssa_total == 0 due to underflow
231              do jg = 1,ng
232                if (od_total(jg) > 0.0_jprb) then
233                  scat_od = ssa(jg,jlev,jcol)*od(jg,jlev,jcol) &
234                       &     + ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
235                       &     *  od_cloud_new(jg)
236                  ssa_total(jg) = scat_od / od_total(jg)
237                  if (scat_od > 0.0_jprb) then
238                    g_total(jg) = (g(jg,jlev,jcol)*ssa(jg,jlev,jcol)*od(jg,jlev,jcol) &
239                         &     +   g_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
240                         &     * ssa_cloud(config%i_band_from_reordered_g_sw(jg),jlev,jcol) &
241                         &     *  od_cloud_new(jg)) &
242                         &     / scat_od
243                  end if
244                end if
245              end do
246
247              ! Apply delta-Eddington scaling to the cloud-aerosol-gas
248              ! mixture
249              if (config%do_sw_delta_scaling_with_gases) then
250                call delta_eddington(od_total, ssa_total, g_total)
251              end if
252
253             ! Compute cloudy-sky reflectance, transmittance etc at
254              ! each model level
255              call calc_two_stream_gammas_sw(ng, &
256                   &  cos_sza, ssa_total, g_total, &
257                   &  gamma1, gamma2, gamma3)
258
259              call calc_reflectance_transmittance_sw(ng, &
260                   &  cos_sza, od_total, ssa_total, &
261                   &  gamma1, gamma2, gamma3, &
262                   &  reflectance(:,jlev), transmittance(:,jlev), &
263                   &  ref_dir(:,jlev), trans_dir_diff(:,jlev), &
264                   &  trans_dir_dir(:,jlev) )
265
266            else
267              ! Clear-sky layer: copy over clear-sky values
268              reflectance(:,jlev) = ref_clear(:,jlev)
269              transmittance(:,jlev) = trans_clear(:,jlev)
270              ref_dir(:,jlev) = ref_dir_clear(:,jlev)
271              trans_dir_diff(:,jlev) = trans_dir_diff_clear(:,jlev)
272              trans_dir_dir(:,jlev) = trans_dir_dir_clear(:,jlev)
273            end if
274          end do
275           
276          ! Use adding method to compute fluxes for an overcast sky
277          call adding_ica_sw(ng, nlev, incoming_sw(:,jcol), &
278               &  albedo_diffuse(:,jcol), albedo_direct(:,jcol), spread(cos_sza,1,ng), &
279               &  reflectance, transmittance, ref_dir, trans_dir_diff, &
280               &  trans_dir_dir, flux_up, flux_dn_diffuse, flux_dn_direct)
281         
282          ! Store overcast broadband fluxes
283          flux%sw_up(jcol,:) = sum(flux_up,1)
284          if (allocated(flux%sw_dn_direct)) then
285            flux%sw_dn_direct(jcol,:) = sum(flux_dn_direct,1)
286            flux%sw_dn(jcol,:) = sum(flux_dn_diffuse,1) &
287                 &  + flux%sw_dn_direct(jcol,:)
288          else
289            flux%sw_dn(jcol,:) = sum(flux_dn_diffuse,1) &
290                 &  + sum(flux_dn_direct,1)
291          end if
292
293          ! Cloudy flux profiles currently assume completely overcast
294          ! skies; perform weighted average with clear-sky profile
295          flux%sw_up(jcol,:) =  total_cloud_cover *flux%sw_up(jcol,:) &
296               &  + (1.0_jprb - total_cloud_cover)*flux%sw_up_clear(jcol,:)
297          flux%sw_dn(jcol,:) =  total_cloud_cover *flux%sw_dn(jcol,:) &
298               &  + (1.0_jprb - total_cloud_cover)*flux%sw_dn_clear(jcol,:)
299          if (allocated(flux%sw_dn_direct)) then
300            flux%sw_dn_direct(jcol,:) = total_cloud_cover *flux%sw_dn_direct(jcol,:) &
301                 &  + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_clear(jcol,:)
302          end if
303          ! Likewise for surface spectral fluxes
304          flux%sw_dn_diffuse_surf_g(:,jcol) = flux_dn_diffuse(:,nlev+1)
305          flux%sw_dn_direct_surf_g(:,jcol)  = flux_dn_direct(:,nlev+1)
306          flux%sw_dn_diffuse_surf_g(:,jcol) = total_cloud_cover *flux%sw_dn_diffuse_surf_g(:,jcol) &
307               &     + (1.0_jprb - total_cloud_cover)*flux%sw_dn_diffuse_surf_clear_g(:,jcol)
308          flux%sw_dn_direct_surf_g(:,jcol) = total_cloud_cover *flux%sw_dn_direct_surf_g(:,jcol) &
309               &     + (1.0_jprb - total_cloud_cover)*flux%sw_dn_direct_surf_clear_g(:,jcol)
310         
311        else
312          ! No cloud in profile and clear-sky fluxes already
313          ! calculated: copy them over
314          flux%sw_up(jcol,:) = flux%sw_up_clear(jcol,:)
315          flux%sw_dn(jcol,:) = flux%sw_dn_clear(jcol,:)
316          if (allocated(flux%sw_dn_direct)) then
317            flux%sw_dn_direct(jcol,:) = flux%sw_dn_direct_clear(jcol,:)
318          end if
319          flux%sw_dn_diffuse_surf_g(:,jcol) = flux%sw_dn_diffuse_surf_clear_g(:,jcol)
320          flux%sw_dn_direct_surf_g(:,jcol)  = flux%sw_dn_direct_surf_clear_g(:,jcol)
321
322        end if ! Cloud is present in profile
323
324      else
325        ! Set fluxes to zero if sun is below the horizon
326        flux%sw_up(jcol,:) = 0.0_jprb
327        flux%sw_dn(jcol,:) = 0.0_jprb
328        if (allocated(flux%sw_dn_direct)) then
329          flux%sw_dn_direct(jcol,:) = 0.0_jprb
330        end if
331        flux%sw_up_clear(jcol,:) = 0.0_jprb
332        flux%sw_dn_clear(jcol,:) = 0.0_jprb
333        if (allocated(flux%sw_dn_direct_clear)) then
334          flux%sw_dn_direct_clear(jcol,:) = 0.0_jprb
335        end if
336        flux%sw_dn_diffuse_surf_g(:,jcol) = 0.0_jprb
337        flux%sw_dn_direct_surf_g(:,jcol)  = 0.0_jprb
338        flux%sw_dn_diffuse_surf_clear_g(:,jcol) = 0.0_jprb
339        flux%sw_dn_direct_surf_clear_g(:,jcol)  = 0.0_jprb
340      end if ! Sun above horizon
341
342    end do ! Loop over columns
343
344    if (lhook) call dr_hook('radiation_mcica_sw:solver_mcica_sw',1,hook_handle)
345   
346  end subroutine solver_mcica_sw
347
348end module radiation_mcica_sw
Note: See TracBrowser for help on using the repository browser.