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

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

Online implementation of the radiative transfer code ECRAD in the LMDZ model.

  • Inclusion of the ecrad directory containing the sources of the ECRAD code
    • interface routine : radiation_scheme.F90
  • Adaptation of compilation scripts :
    • compilation under CPP key CPP_ECRAD
    • compilation with option "-rad ecard" or "-ecard true"
    • The "-rad old/rtm/ecran" build option will need to replace the "-rrtm true" and "-ecrad true" options in the future.
  • Runing LMDZ simulations with ecrad, you need :
    • logical key iflag_rrtm = 2 in physiq.def
    • namelist_ecrad (DefLists?)
    • the directory "data" containing the configuration files is temporarily placed in ../libfphylmd/ecrad/
  • Compilation and execution are tested in the 1D case. The repository under svn would allow to continue the implementation work: tests, verification of the results, ...
File size: 6.1 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
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(jprb) :: 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    do jlev = nlev,1,-1
101      ! Next loop over columns. We could do this by indexing the
102      ! entire inner dimension as follows, e.g. for the first line:
103      !   inv_denominator(:,jlev) = 1.0_jprb / (1.0_jprb-albedo(:,jlev+1)*reflectance(:,jlev))
104      ! and similarly for subsequent lines, but this slows down the
105      ! routine by a factor of 2!  Rather, we do it with an explicit
106      ! loop.
107      do jcol = 1,ncol
108        ! Lacis and Hansen (1974) Eq 33, Shonk & Hogan (2008) Eq 10:
109        inv_denominator(jcol,jlev) = 1.0_jprb / (1.0_jprb-albedo(jcol,jlev+1)*reflectance(jcol,jlev))
110        ! Shonk & Hogan (2008) Eq 9, Petty (2006) Eq 13.81:
111        albedo(jcol,jlev) = reflectance(jcol,jlev) + transmittance(jcol,jlev) * transmittance(jcol,jlev) &
112             &                                     * albedo(jcol,jlev+1) * inv_denominator(jcol,jlev)
113        ! Shonk & Hogan (2008) Eq 11:
114        source(jcol,jlev) = ref_dir(jcol,jlev)*flux_dn_direct(jcol,jlev) &
115             &  + transmittance(jcol,jlev)*(source(jcol,jlev+1) &
116             &        + albedo(jcol,jlev+1)*trans_dir_diff(jcol,jlev)*flux_dn_direct(jcol,jlev)) &
117             &  * inv_denominator(jcol,jlev)
118      end do
119    end do
120
121    ! At top-of-atmosphere there is no diffuse downwelling radiation
122    flux_dn_diffuse(:,1) = 0.0_jprb
123
124    ! At top-of-atmosphere, all upwelling radiation is due to
125    ! scattering by the direct beam below that level
126    flux_up(:,1) = source(:,1)
127
128    ! Work back down through the atmosphere computing the fluxes at
129    ! each half-level
130    do jlev = 1,nlev
131      do jcol = 1,ncol
132        ! Shonk & Hogan (2008) Eq 14 (after simplification):
133        flux_dn_diffuse(jcol,jlev+1) &
134             &  = (transmittance(jcol,jlev)*flux_dn_diffuse(jcol,jlev) &
135             &     + reflectance(jcol,jlev)*source(jcol,jlev+1) &
136             &     + trans_dir_diff(jcol,jlev)*flux_dn_direct(jcol,jlev)) * inv_denominator(jcol,jlev)
137        ! Shonk & Hogan (2008) Eq 12:
138        flux_up(jcol,jlev+1) = albedo(jcol,jlev+1)*flux_dn_diffuse(jcol,jlev+1) &
139             &            + source(jcol,jlev+1)
140        flux_dn_direct(jcol,jlev) = flux_dn_direct(jcol,jlev)*cos_sza(jcol)
141      end do
142    end do
143    flux_dn_direct(:,nlev+1) = flux_dn_direct(:,nlev+1)*cos_sza
144
145    if (lhook) call dr_hook('radiation_adding_ica_sw:adding_ica_sw',1,hook_handle)
146
147  end subroutine adding_ica_sw
148
149end module radiation_adding_ica_sw
Note: See TracBrowser for help on using the repository browser.