source: LMDZ6/trunk/libf/phylmd/ecrad/radiation/radiation_adding_ica_sw.F90 @ 5413

Last change on this file since 5413 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.2 KB
Line 
1! radiation_adding_ica_sw.F90 - Shortwave adding method in independent column approximation
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-10-23  R. Hogan  Renamed single-character variables
17
18module radiation_adding_ica_sw
19
20  public
21
22contains
23
24  subroutine adding_ica_sw(ncol, nlev, incoming_toa, &
25       &  albedo_surf_diffuse, albedo_surf_direct, cos_sza, &
26       &  reflectance, transmittance, ref_dir, trans_dir_diff, trans_dir_dir, &
27       &  flux_up, flux_dn_diffuse, flux_dn_direct)
28
29    use parkind1, only           : jprb
30    use yomhook,  only           : lhook, dr_hook, jphook
31
32    implicit none
33
34    ! Inputs
35    integer, intent(in) :: ncol ! number of columns (may be spectral intervals)
36    integer, intent(in) :: nlev ! number of levels
37
38    ! Incoming downwelling solar radiation at top-of-atmosphere (W m-2)
39    real(jprb), intent(in),  dimension(ncol)         :: incoming_toa
40
41    ! Surface albedo to diffuse and direct radiation
42    real(jprb), intent(in),  dimension(ncol)         :: albedo_surf_diffuse, &
43         &                                              albedo_surf_direct
44
45    ! Cosine of the solar zenith angle
46    real(jprb), intent(in),  dimension(ncol)         :: cos_sza
47
48    ! Diffuse reflectance and transmittance of each layer
49    real(jprb), intent(in),  dimension(ncol, nlev)   :: reflectance, transmittance
50
51    ! Fraction of direct-beam solar radiation entering the top of a
52    ! layer that is reflected back up or scattered forward into the
53    ! diffuse stream at the base of the layer
54    real(jprb), intent(in),  dimension(ncol, nlev)   :: ref_dir, trans_dir_diff
55
56    ! Direct transmittance, i.e. fraction of direct beam that
57    ! penetrates a layer without being scattered or absorbed
58    real(jprb), intent(in),  dimension(ncol, nlev)   :: trans_dir_dir
59
60    ! Resulting fluxes (W m-2) at half-levels: diffuse upwelling,
61    ! diffuse downwelling and direct downwelling
62    real(jprb), intent(out), dimension(ncol, nlev+1) :: flux_up, flux_dn_diffuse, &
63         &                                              flux_dn_direct
64   
65    ! Albedo of the entire earth/atmosphere system below each half
66    ! level
67    real(jprb), dimension(ncol, nlev+1) :: albedo
68
69    ! Upwelling radiation at each half-level due to scattering of the
70    ! direct beam below that half-level (W m-2)
71    real(jprb), dimension(ncol, nlev+1) :: source
72
73    ! Equal to 1/(1-albedo*reflectance)
74    real(jprb), dimension(ncol, nlev)   :: inv_denominator
75
76    ! Loop index for model level and column
77    integer :: jlev, jcol
78
79    real(jphook) :: hook_handle
80
81    if (lhook) call dr_hook('radiation_adding_ica_sw:adding_ica_sw',0,hook_handle)
82
83    ! Compute profile of direct (unscattered) solar fluxes at each
84    ! half-level by working down through the atmosphere
85    flux_dn_direct(:,1) = incoming_toa
86    do jlev = 1,nlev
87      flux_dn_direct(:,jlev+1) = flux_dn_direct(:,jlev)*trans_dir_dir(:,jlev)
88    end do
89
90    albedo(:,nlev+1) = albedo_surf_diffuse
91
92    ! At the surface, the direct solar beam is reflected back into the
93    ! diffuse stream
94    source(:,nlev+1) = albedo_surf_direct * flux_dn_direct(:,nlev+1) * cos_sza
95
96    ! Work back up through the atmosphere and compute the albedo of
97    ! the entire earth/atmosphere system below that half-level, and
98    ! also the "source", which is the upwelling flux due to direct
99    ! radiation that is scattered below that level
100! Added for DWD (2020)
101!NEC$ outerloop_unroll(8)
102    do jlev = nlev,1,-1
103      ! Next loop over columns. We could do this by indexing the
104      ! entire inner dimension as follows, e.g. for the first line:
105      !   inv_denominator(:,jlev) = 1.0_jprb / (1.0_jprb-albedo(:,jlev+1)*reflectance(:,jlev))
106      ! and similarly for subsequent lines, but this slows down the
107      ! routine by a factor of 2!  Rather, we do it with an explicit
108      ! loop.
109      do jcol = 1,ncol
110        ! Lacis and Hansen (1974) Eq 33, Shonk & Hogan (2008) Eq 10:
111        inv_denominator(jcol,jlev) = 1.0_jprb / (1.0_jprb-albedo(jcol,jlev+1)*reflectance(jcol,jlev))
112        ! Shonk & Hogan (2008) Eq 9, Petty (2006) Eq 13.81:
113        albedo(jcol,jlev) = reflectance(jcol,jlev) + transmittance(jcol,jlev) * transmittance(jcol,jlev) &
114             &                                     * albedo(jcol,jlev+1) * inv_denominator(jcol,jlev)
115        ! Shonk & Hogan (2008) Eq 11:
116        source(jcol,jlev) = ref_dir(jcol,jlev)*flux_dn_direct(jcol,jlev) &
117             &  + transmittance(jcol,jlev)*(source(jcol,jlev+1) &
118             &        + albedo(jcol,jlev+1)*trans_dir_diff(jcol,jlev)*flux_dn_direct(jcol,jlev)) &
119             &  * inv_denominator(jcol,jlev)
120      end do
121    end do
122
123    ! At top-of-atmosphere there is no diffuse downwelling radiation
124    flux_dn_diffuse(:,1) = 0.0_jprb
125
126    ! At top-of-atmosphere, all upwelling radiation is due to
127    ! scattering by the direct beam below that level
128    flux_up(:,1) = source(:,1)
129
130    ! Work back down through the atmosphere computing the fluxes at
131    ! each half-level
132! Added for DWD (2020)
133!NEC$ outerloop_unroll(8)
134    do jlev = 1,nlev
135      do jcol = 1,ncol
136        ! Shonk & Hogan (2008) Eq 14 (after simplification):
137        flux_dn_diffuse(jcol,jlev+1) &
138             &  = (transmittance(jcol,jlev)*flux_dn_diffuse(jcol,jlev) &
139             &     + reflectance(jcol,jlev)*source(jcol,jlev+1) &
140             &     + trans_dir_diff(jcol,jlev)*flux_dn_direct(jcol,jlev)) * inv_denominator(jcol,jlev)
141        ! Shonk & Hogan (2008) Eq 12:
142        flux_up(jcol,jlev+1) = albedo(jcol,jlev+1)*flux_dn_diffuse(jcol,jlev+1) &
143             &            + source(jcol,jlev+1)
144        flux_dn_direct(jcol,jlev) = flux_dn_direct(jcol,jlev)*cos_sza(jcol)
145      end do
146    end do
147    flux_dn_direct(:,nlev+1) = flux_dn_direct(:,nlev+1)*cos_sza
148
149    if (lhook) call dr_hook('radiation_adding_ica_sw:adding_ica_sw',1,hook_handle)
150
151  end subroutine adding_ica_sw
152
153end module radiation_adding_ica_sw
Note: See TracBrowser for help on using the repository browser.