source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_cloudless_lw.F90 @ 5418

Last change on this file since 5418 was 4773, checked in by idelkadi, 12 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: 6.7 KB
Line 
1! radiation_cloudless_lw.F90 - Longwave homogeneous cloudless solver
2!
3! (C) Copyright 2019- 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
16module radiation_cloudless_lw
17
18public :: solver_cloudless_lw
19
20contains
21
22  !---------------------------------------------------------------------
23  ! Longwave homogeneous solver containing no clouds
24  subroutine solver_cloudless_lw(nlev,istartcol,iendcol, &
25       &  config, od, ssa, g, planck_hl, emission, albedo, flux)
26
27    use parkind1, only           : jprb
28    use yomhook,  only           : lhook, dr_hook, jphook
29
30    use radiation_config, only         : config_type
31    use radiation_flux, only           : flux_type, indexed_sum_profile
32    use radiation_two_stream, only     : calc_two_stream_gammas_lw, &
33         &                               calc_reflectance_transmittance_lw, &
34         &                               calc_no_scattering_transmittance_lw
35    use radiation_adding_ica_lw, only  : adding_ica_lw, calc_fluxes_no_scattering_lw
36    use radiation_lw_derivatives, only : calc_lw_derivatives_ica
37
38    implicit none
39
40    ! Inputs
41    integer, intent(in) :: nlev               ! number of model levels
42    integer, intent(in) :: istartcol, iendcol ! range of columns to process
43    type(config_type),        intent(in) :: config
44
45    ! Gas and aerosol optical depth, single-scattering albedo and
46    ! asymmetry factor at each longwave g-point
47    real(jprb), intent(in), dimension(config%n_g_lw, nlev, istartcol:iendcol) :: &
48         &  od
49    real(jprb), intent(in), dimension(config%n_g_lw_if_scattering, nlev, istartcol:iendcol) :: &
50         &  ssa, g
51
52    ! Planck function at each half-level and the surface
53    real(jprb), intent(in), dimension(config%n_g_lw,nlev+1,istartcol:iendcol) :: &
54         &  planck_hl
55 
56    ! Emission (Planck*emissivity) and albedo (1-emissivity) at the
57    ! surface at each longwave g-point
58    real(jprb), intent(in), dimension(config%n_g_lw, istartcol:iendcol) &
59         &  :: emission, albedo
60
61    ! Output
62    type(flux_type), intent(inout):: flux
63
64    ! Local variables
65
66    ! Diffuse reflectance and transmittance for each layer in clear
67    ! and all skies
68    real(jprb), dimension(config%n_g_lw, nlev) :: reflectance, transmittance
69
70    ! Emission by a layer into the upwelling or downwelling diffuse
71    ! streams, in clear and all skies
72    real(jprb), dimension(config%n_g_lw, nlev) :: source_up, source_dn
73
74    ! Fluxes per g point
75    real(jprb), dimension(config%n_g_lw, nlev+1) :: flux_up, flux_dn
76
77    ! Combined optical depth, single scattering albedo and asymmetry
78    ! factor
79    real(jprb), dimension(config%n_g_lw) :: ssa_total, g_total
80
81    ! Two-stream coefficients
82    real(jprb), dimension(config%n_g_lw) :: gamma1, gamma2
83
84    ! Number of g points
85    integer :: ng
86
87    ! Loop indices for level and column
88    integer :: jlev, jcol
89
90    real(jphook) :: hook_handle
91
92    if (lhook) call dr_hook('radiation_cloudless_lw:solver_cloudless_lw',0,hook_handle)
93
94    ng = config%n_g_lw
95
96    ! Loop through columns
97    do jcol = istartcol,iendcol
98
99      ! Compute the reflectance and transmittance of all layers,
100      ! neglecting clouds
101      do jlev = 1,nlev
102        if (config%do_lw_aerosol_scattering) then
103          ! Scattering case: first compute clear-sky reflectance,
104          ! transmittance etc at each model level
105          ssa_total = ssa(:,jlev,jcol)
106          g_total   = g(:,jlev,jcol)
107          call calc_two_stream_gammas_lw(ng, ssa_total, g_total, &
108               &  gamma1, gamma2)
109          call calc_reflectance_transmittance_lw(ng, &
110               &  od(:,jlev,jcol), gamma1, gamma2, &
111               &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1,jcol), &
112               &  reflectance(:,jlev), transmittance(:,jlev), &
113               &  source_up(:,jlev), source_dn(:,jlev))
114        else
115          ! Non-scattering case: use simpler functions for
116          ! transmission and emission
117          call calc_no_scattering_transmittance_lw(ng, od(:,jlev,jcol), &
118               &  planck_hl(:,jlev,jcol), planck_hl(:,jlev+1, jcol), &
119               &  transmittance(:,jlev), source_up(:,jlev), source_dn(:,jlev))         
120          ! Ensure that clear-sky reflectance is zero
121          reflectance(:,jlev) = 0.0_jprb
122        end if
123      end do
124
125      if (config%do_lw_aerosol_scattering) then
126        ! Then use adding method to compute fluxes
127        call adding_ica_lw(ng, nlev, &
128             &  reflectance, transmittance, source_up, source_dn, &
129             &  emission(:,jcol), albedo(:,jcol), &
130             &  flux_up, flux_dn)
131      else
132        ! Simpler down-then-up method to compute fluxes
133        call calc_fluxes_no_scattering_lw(ng, nlev, &
134             &  transmittance, source_up, source_dn, &
135             &  emission(:,jcol), albedo(:,jcol), &
136             &  flux_up, flux_dn)
137         
138      end if
139
140      ! Sum over g-points to compute broadband fluxes
141      flux%lw_up(jcol,:) = sum(flux_up,1)
142      flux%lw_dn(jcol,:) = sum(flux_dn,1)
143      ! Store surface spectral downwelling fluxes
144      flux%lw_dn_surf_g(:,jcol) = flux_dn(:,nlev+1)
145
146      ! Save the spectral fluxes if required
147      if (config%do_save_spectral_flux) then
148        call indexed_sum_profile(flux_up, config%i_spec_from_reordered_g_lw, &
149             &                   flux%lw_up_band(:,jcol,:))
150        call indexed_sum_profile(flux_dn, config%i_spec_from_reordered_g_lw, &
151             &                   flux%lw_dn_band(:,jcol,:))
152      end if
153
154      if (config%do_clear) then
155        ! Clear-sky calculations are equal to all-sky for this solver:
156        ! copy fluxes over
157        flux%lw_up_clear(jcol,:) = flux%lw_up(jcol,:)
158        flux%lw_dn_clear(jcol,:) = flux%lw_dn(jcol,:)
159        flux%lw_dn_surf_clear_g(:,jcol) = flux%lw_dn_surf_g(:,jcol)
160        if (config%do_save_spectral_flux) then
161          flux%lw_up_clear_band(:,jcol,:) = flux%lw_up_band(:,jcol,:)
162          flux%lw_dn_clear_band(:,jcol,:) = flux%lw_dn_band(:,jcol,:)
163        end if
164      end if
165
166      ! Compute the longwave derivatives needed by Hogan and Bozzo
167      ! (2015) approximate radiation update scheme
168      if (config%do_lw_derivatives) then
169        call calc_lw_derivatives_ica(ng, nlev, jcol, transmittance, flux_up(:,nlev+1), &
170             &                       flux%lw_derivatives)
171       end if
172
173    end do
174
175    if (lhook) call dr_hook('radiation_cloudless_lw:solver_cloudless_lw',1,hook_handle)
176   
177  end subroutine solver_cloudless_lw
178
179end module radiation_cloudless_lw
Note: See TracBrowser for help on using the repository browser.