source: LMDZ6/branches/cirrus/libf/phylmd/ecrad/radiation/radiation_mcica_lw.F90

Last change on this file was 5202, checked in by Laurent Fairhead, 8 weeks ago

Updating cirrus branch to trunk revision 5171

File size: 15.2 KB
Line 
1! radiation_mcica_lw.F90 - Monte-Carlo Independent Column Approximation longtwave 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 emission/albedo rather than planck/emissivity
17!   2017-04-22  R. Hogan  Store surface fluxes at all g-points
18!   2017-07-12  R. Hogan  Call fast adding method if only clouds scatter
19!   2017-10-23  R. Hogan  Renamed single-character variables
20
21module radiation_mcica_lw
22
23  public
24
25contains
26
27  !---------------------------------------------------------------------
28  ! Longwave Monte Carlo Independent Column Approximation
29  ! (McICA). This implementation performs a clear-sky and a cloudy-sky
30  ! calculation, and then weights the two to get the all-sky fluxes
31  ! according to the total cloud cover. This method reduces noise for
32  ! low cloud cover situations, and exploits the clear-sky
33  ! calculations that are usually performed for diagnostic purposes
34  ! simultaneously. The cloud generator has been carefully written
35  ! such that the stochastic cloud field satisfies the prescribed
36  ! overlap parameter accounting for this weighting.
37  subroutine solver_mcica_lw(nlev,istartcol,iendcol, &
38       &  config, single_level, cloud, &
39       &  od, ssa, g, od_cloud, ssa_cloud, g_cloud, planck_hl, &
40       &  emission, albedo, &
41       &  flux)
42
43    use parkind1, only           : jprb
44    use yomhook,  only           : lhook, dr_hook, jphook
45
46    use radiation_io,   only           : nulerr, radiation_abort
47    use radiation_config, only         : config_type
48    use radiation_single_level, only   : single_level_type
49    use radiation_cloud, only          : cloud_type
50    use radiation_flux, only           : flux_type
51    use radiation_two_stream, only     : calc_ref_trans_lw, &
52         &                               calc_no_scattering_transmittance_lw
53    use radiation_adding_ica_lw, only  : adding_ica_lw, fast_adding_ica_lw, &
54         &                               calc_fluxes_no_scattering_lw
55    use radiation_lw_derivatives, only : calc_lw_derivatives_ica, modify_lw_derivatives_ica
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 longwave g-point
69    real(jprb), intent(in), dimension(config%n_g_lw, nlev, istartcol:iendcol) :: &
70         &  od
71    real(jprb), intent(in), dimension(config%n_g_lw_if_scattering, nlev, istartcol:iendcol) :: &
72         &  ssa, g
73
74    ! Cloud and precipitation optical depth, single-scattering albedo and
75    ! asymmetry factor in each longwave band
76    real(jprb), intent(in), dimension(config%n_bands_lw,nlev,istartcol:iendcol)   :: &
77         &  od_cloud
78    real(jprb), intent(in), dimension(config%n_bands_lw_if_scattering, &
79         &  nlev,istartcol:iendcol) :: ssa_cloud, g_cloud
80
81    ! Planck function at each half-level and the surface
82    real(jprb), intent(in), dimension(config%n_g_lw,nlev+1,istartcol:iendcol) :: &
83         &  planck_hl
84
85    ! Emission (Planck*emissivity) and albedo (1-emissivity) at the
86    ! surface at each longwave g-point
87    real(jprb), intent(in), dimension(config%n_g_lw, istartcol:iendcol) :: emission, albedo
88
89    ! Output
90    type(flux_type), intent(inout):: flux
91
92    ! Local variables
93
94    ! Diffuse reflectance and transmittance for each layer in clear
95    ! and all skies
96    real(jprb), dimension(config%n_g_lw, nlev) :: ref_clear, trans_clear, reflectance, transmittance
97
98    ! Emission by a layer into the upwelling or downwelling diffuse
99    ! streams, in clear and all skies
100    real(jprb), dimension(config%n_g_lw, nlev) :: source_up_clear, source_dn_clear, source_up, source_dn
101
102    ! Fluxes per g point
103    real(jprb), dimension(config%n_g_lw, nlev+1) :: flux_up, flux_dn
104    real(jprb), dimension(config%n_g_lw, nlev+1) :: flux_up_clear, flux_dn_clear
105
106    ! Combined gas+aerosol+cloud optical depth, single scattering
107    ! albedo and asymmetry factor
108    real(jprb), dimension(config%n_g_lw) :: od_total, ssa_total, g_total
109
110    ! Combined scattering optical depth
111    real(jprb) :: scat_od, scat_od_total(config%n_g_lw)
112
113    ! Optical depth scaling from the cloud generator, zero indicating
114    ! clear skies
115    real(jprb), dimension(config%n_g_lw,nlev) :: od_scaling
116
117    ! Modified optical depth after McICA scaling to represent cloud
118    ! inhomogeneity
119    real(jprb), dimension(config%n_g_lw) :: od_cloud_new
120
121    ! Total cloud cover output from the cloud generator
122    real(jprb) :: total_cloud_cover
123
124    ! Identify clear-sky layers
125    logical :: is_clear_sky_layer(nlev)
126
127    ! Index of the highest cloudy layer
128    integer :: i_cloud_top
129
130    ! Number of g points
131    integer :: ng
132
133    ! Loop indices for level, column and g point
134    integer :: jlev, jcol, jg
135
136    real(jphook) :: hook_handle
137
138    if (lhook) call dr_hook('radiation_mcica_lw:solver_mcica_lw',0,hook_handle)
139
140    if (.not. config%do_clear) then
141      write(nulerr,'(a)') '*** Error: longwave McICA requires clear-sky calculation to be performed'
142      call radiation_abort()     
143    end if
144
145    ng = config%n_g_lw
146
147    ! Loop through columns
148    do jcol = istartcol,iendcol
149
150      ! Clear-sky calculation
151      if (config%do_lw_aerosol_scattering) then
152        ! Scattering case: first compute clear-sky reflectance,
153        ! transmittance etc at each model level
154        call calc_ref_trans_lw(ng*nlev, &
155             &  od(:,:,jcol), ssa(:,:,jcol), g(:,:,jcol), &
156             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1,jcol), &
157             &  ref_clear, trans_clear, &
158             &  source_up_clear, source_dn_clear)
159        ! Then use adding method to compute fluxes
160        call adding_ica_lw(ng, nlev, &
161             &  ref_clear, trans_clear, source_up_clear, source_dn_clear, &
162             &  emission(:,jcol), albedo(:,jcol), &
163             &  flux_up_clear, flux_dn_clear)
164      else
165        ! Non-scattering case: use simpler functions for
166        ! transmission and emission
167        call calc_no_scattering_transmittance_lw(ng*nlev, od(:,:,jcol), &
168             &  planck_hl(:,1:nlev,jcol), planck_hl(:,2:nlev+1, jcol), &
169             &  trans_clear, source_up_clear, source_dn_clear)
170        ! Ensure that clear-sky reflectance is zero since it may be
171        ! used in cloudy-sky case
172        ref_clear = 0.0_jprb
173        ! Simpler down-then-up method to compute fluxes
174        call calc_fluxes_no_scattering_lw(ng, nlev, &
175             &  trans_clear, source_up_clear, source_dn_clear, &
176             &  emission(:,jcol), albedo(:,jcol), &
177             &  flux_up_clear, flux_dn_clear)       
178      end if
179
180      ! Sum over g-points to compute broadband fluxes
181      flux%lw_up_clear(jcol,:) = sum(flux_up_clear,1)
182      flux%lw_dn_clear(jcol,:) = sum(flux_dn_clear,1)
183      ! Store surface spectral downwelling fluxes
184      flux%lw_dn_surf_clear_g(:,jcol) = flux_dn_clear(:,nlev+1)
185
186      ! Do cloudy-sky calculation; add a prime number to the seed in
187      ! the longwave
188      call cloud_generator(ng, nlev, config%i_overlap_scheme, &
189           &  single_level%iseed(jcol) + 997, &
190           &  config%cloud_fraction_threshold, &
191           &  cloud%fraction(jcol,:), cloud%overlap_param(jcol,:), &
192           &  config%cloud_inhom_decorr_scaling, cloud%fractional_std(jcol,:), &
193           &  config%pdf_sampler, od_scaling, total_cloud_cover, &
194           &  use_beta_overlap=config%use_beta_overlap, &
195           &  use_vectorizable_generator=config%use_vectorizable_generator)
196     
197      ! Store total cloud cover
198      flux%cloud_cover_lw(jcol) = total_cloud_cover
199     
200      if (total_cloud_cover >= config%cloud_fraction_threshold) then
201        ! Total-sky calculation
202
203        is_clear_sky_layer = .true.
204        i_cloud_top = nlev+1
205        do jlev = 1,nlev
206          ! Compute combined gas+aerosol+cloud optical properties
207          if (cloud%fraction(jcol,jlev) >= config%cloud_fraction_threshold) then
208            is_clear_sky_layer(jlev) = .false.
209            ! Get index to the first cloudy layer from the top
210            if (i_cloud_top > jlev) then
211              i_cloud_top = jlev
212            end if
213
214            do jg = 1,ng
215              od_cloud_new(jg) = od_scaling(jg,jlev) &
216                 &  * od_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol)
217              od_total(jg)  = od(jg,jlev,jcol) + od_cloud_new(jg)
218              ssa_total(jg) = 0.0_jprb
219              g_total(jg)   = 0.0_jprb
220            end do
221
222            if (config%do_lw_cloud_scattering) then
223              ! Scattering case: calculate reflectance and
224              ! transmittance at each model level
225
226              if (config%do_lw_aerosol_scattering) then
227                ! In single precision we need to protect against the
228                ! case that od_total > 0.0 and ssa_total > 0.0 but
229                ! od_total*ssa_total == 0 due to underflow
230                do jg = 1,ng
231                  if (od_total(jg) > 0.0_jprb) then
232                    scat_od_total(jg) = ssa(jg,jlev,jcol)*od(jg,jlev,jcol) &
233                     &     + ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) &
234                     &     *  od_cloud_new(jg)
235                    ssa_total(jg) = scat_od_total(jg) / od_total(jg)
236
237                    if (scat_od_total(jg) > 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_lw(jg),jlev,jcol) &
240                         &     * ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) &
241                         &     *  od_cloud_new(jg)) &
242                         &     / scat_od_total(jg)
243                    end if
244                  end if
245                end do
246
247              else
248
249                do jg = 1,ng
250                  if (od_total(jg) > 0.0_jprb) then
251                    scat_od = ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) &
252                         &     * od_cloud_new(jg)
253                    ssa_total(jg) = scat_od / od_total(jg)
254                    if (scat_od > 0.0_jprb) then
255                      g_total(jg) = g_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) &
256                           &     * ssa_cloud(config%i_band_from_reordered_g_lw(jg),jlev,jcol) &
257                           &     *  od_cloud_new(jg) / scat_od
258                    end if
259                  end if
260                end do
261
262              end if
263           
264              ! Compute cloudy-sky reflectance, transmittance etc at
265              ! each model level
266              call calc_ref_trans_lw(ng, &
267                   &  od_total, ssa_total, g_total, &
268                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1,jcol), &
269                   &  reflectance(:,jlev), transmittance(:,jlev), &
270                   &  source_up(:,jlev), source_dn(:,jlev))
271            else
272              ! No-scattering case: use simpler functions for
273              ! transmission and emission
274              call calc_no_scattering_transmittance_lw(ng, od_total, &
275                   &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1, jcol), &
276                   &  transmittance(:,jlev), source_up(:,jlev), source_dn(:,jlev))
277            end if
278
279          else
280            ! Clear-sky layer: copy over clear-sky values
281            reflectance(:,jlev) = ref_clear(:,jlev)
282            transmittance(:,jlev) = trans_clear(:,jlev)
283            source_up(:,jlev) = source_up_clear(:,jlev)
284            source_dn(:,jlev) = source_dn_clear(:,jlev)
285          end if
286        end do
287       
288        if (config%do_lw_aerosol_scattering) then
289          ! Use adding method to compute fluxes for an overcast sky,
290          ! allowing for scattering in all layers
291          call adding_ica_lw(ng, nlev, reflectance, transmittance, source_up, source_dn, &
292               &  emission(:,jcol), albedo(:,jcol), &
293               &  flux_up, flux_dn)
294        else if (config%do_lw_cloud_scattering) then
295          ! Use adding method to compute fluxes but optimize for the
296          ! presence of clear-sky layers
297          call fast_adding_ica_lw(ng, nlev, reflectance, transmittance, source_up, source_dn, &
298               &  emission(:,jcol), albedo(:,jcol), &
299               &  is_clear_sky_layer, i_cloud_top, flux_dn_clear, &
300               &  flux_up, flux_dn)
301        else
302          ! Simpler down-then-up method to compute fluxes
303          call calc_fluxes_no_scattering_lw(ng, nlev, &
304               &  transmittance, source_up, source_dn, emission(:,jcol), albedo(:,jcol), &
305               &  flux_up, flux_dn)
306        end if
307       
308        ! Store overcast broadband fluxes
309        flux%lw_up(jcol,:) = sum(flux_up,1)
310        flux%lw_dn(jcol,:) = sum(flux_dn,1)
311
312        ! Cloudy flux profiles currently assume completely overcast
313        ! skies; perform weighted average with clear-sky profile
314        flux%lw_up(jcol,:) =  total_cloud_cover *flux%lw_up(jcol,:) &
315             &  + (1.0_jprb - total_cloud_cover)*flux%lw_up_clear(jcol,:)
316        flux%lw_dn(jcol,:) =  total_cloud_cover *flux%lw_dn(jcol,:) &
317             &  + (1.0_jprb - total_cloud_cover)*flux%lw_dn_clear(jcol,:)
318        ! Store surface spectral downwelling fluxes
319        flux%lw_dn_surf_g(:,jcol) = total_cloud_cover*flux_dn(:,nlev+1) &
320             &  + (1.0_jprb - total_cloud_cover)*flux%lw_dn_surf_clear_g(:,jcol)
321
322        ! Compute the longwave derivatives needed by Hogan and Bozzo
323        ! (2015) approximate radiation update scheme
324        if (config%do_lw_derivatives) then
325          call calc_lw_derivatives_ica(ng, nlev, jcol, transmittance, flux_up(:,nlev+1), &
326               &                       flux%lw_derivatives)
327          if (total_cloud_cover < 1.0_jprb - config%cloud_fraction_threshold) then
328            ! Modify the existing derivative with the contribution from the clear sky
329            call modify_lw_derivatives_ica(ng, nlev, jcol, trans_clear, flux_up_clear(:,nlev+1), &
330                 &                         1.0_jprb-total_cloud_cover, flux%lw_derivatives)
331          end if
332        end if
333
334      else
335        ! No cloud in profile and clear-sky fluxes already
336        ! calculated: copy them over
337        flux%lw_up(jcol,:) = flux%lw_up_clear(jcol,:)
338        flux%lw_dn(jcol,:) = flux%lw_dn_clear(jcol,:)
339        flux%lw_dn_surf_g(:,jcol) = flux%lw_dn_surf_clear_g(:,jcol)
340        if (config%do_lw_derivatives) then
341          call calc_lw_derivatives_ica(ng, nlev, jcol, trans_clear, flux_up_clear(:,nlev+1), &
342               &                       flux%lw_derivatives)
343 
344        end if
345      end if ! Cloud is present in profile
346    end do
347
348    if (lhook) call dr_hook('radiation_mcica_lw:solver_mcica_lw',1,hook_handle)
349   
350  end subroutine solver_mcica_lw
351
352end module radiation_mcica_lw
Note: See TracBrowser for help on using the repository browser.