source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_lw_derivatives.F90 @ 5423

Last change on this file since 5423 was 3908, checked in by idelkadi, 4 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: 9.7 KB
Line 
1! radiation_lw_derivatives.F90 - Compute longwave derivatives for Hogan and Bozzo (2015) method
2!
3! (C) Copyright 2016- 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! This module provides routines to compute the rate of change of
16! broadband upwelling longwave flux at each half level with respect to
17! the surface broadband upwelling flux.  This is done from the surface
18! spectral fluxes and the spectral transmittance of each atmospheric
19! layer, assuming no longwave scattering. The result may be used to
20! perform approximate updates to the longwave flux profile in between
21! calls to the full radiation scheme, accounting for the change in
22! skin temperature, following the method of Hogan and Bozzo (JAMES
23! 2015).  Separate routines are provided for each solver.
24!
25! Note that currently a more approximate calculation is performed from
26! the exact one in Hogan and Bozzo (2015); here we assume that a
27! change in temperature increases the spectral fluxes in proportion,
28! when in reality there is a change in shape of the Planck function in
29! addition to an overall increase in the total emission.
30!
31! Modifications
32!   2017-10-23  R. Hogan  Renamed single-character variables
33
34module radiation_lw_derivatives
35
36  public
37
38contains
39
40  !---------------------------------------------------------------------
41  ! Calculation for the Independent Column Approximation
42  subroutine calc_lw_derivatives_ica(ng, nlev, icol, transmittance, flux_up_surf, lw_derivatives)
43
44    use parkind1, only           : jprb
45    use yomhook,  only           : lhook, dr_hook
46
47    implicit none
48
49    ! Inputs
50    integer,    intent(in) :: ng   ! number of spectral intervals
51    integer,    intent(in) :: nlev ! number of levels
52    integer,    intent(in) :: icol ! Index of column for output
53    real(jprb), intent(in) :: transmittance(ng,nlev)
54    real(jprb), intent(in) :: flux_up_surf(ng) ! Upwelling surface spectral flux (W m-2)
55   
56    ! Output
57    real(jprb), intent(out) :: lw_derivatives(:,:) ! dimensioned (ncol,nlev+1)
58
59    ! Rate of change of spectral flux at a given height with respect
60    ! to the surface value
61    real(jprb) :: lw_derivatives_g(ng)
62
63    integer    :: jlev
64
65    real(jprb) :: hook_handle
66
67    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_ica',0,hook_handle)
68
69    ! Initialize the derivatives at the surface
70    lw_derivatives_g = flux_up_surf / sum(flux_up_surf)
71    lw_derivatives(icol, nlev+1) = 1.0_jprb
72
73    ! Move up through the atmosphere computing the derivatives at each
74    ! half-level
75    do jlev = nlev,1,-1
76      lw_derivatives_g = lw_derivatives_g * transmittance(:,jlev)
77      lw_derivatives(icol,jlev) = sum(lw_derivatives_g)
78    end do
79
80    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_ica',1,hook_handle)
81
82  end subroutine calc_lw_derivatives_ica
83
84
85  !---------------------------------------------------------------------
86  ! Calculation for the Independent Column Approximation
87  subroutine modify_lw_derivatives_ica(ng, nlev, icol, transmittance, &
88       &                               flux_up_surf, weight, lw_derivatives)
89
90    use parkind1, only           : jprb
91    use yomhook,  only           : lhook, dr_hook
92
93    implicit none
94
95    ! Inputs
96    integer,    intent(in) :: ng   ! number of spectral intervals
97    integer,    intent(in) :: nlev ! number of levels
98    integer,    intent(in) :: icol ! Index of column for output
99    real(jprb), intent(in) :: transmittance(ng,nlev)
100    real(jprb), intent(in) :: flux_up_surf(ng) ! Upwelling surface spectral flux (W m-2)
101    real(jprb), intent(in) :: weight ! Weight new values against existing
102   
103    ! Output
104    real(jprb), intent(inout) :: lw_derivatives(:,:) ! dimensioned (ncol,nlev+1)
105
106    ! Rate of change of spectral flux at a given height with respect
107    ! to the surface value
108    real(jprb) :: lw_derivatives_g(ng)
109
110    integer    :: jlev
111
112    real(jprb) :: hook_handle
113
114    if (lhook) call dr_hook('radiation_lw_derivatives:modify_lw_derivatives_ica',0,hook_handle)
115
116    ! Initialize the derivatives at the surface
117    lw_derivatives_g = flux_up_surf / sum(flux_up_surf)
118    ! This value must be 1 so no weighting applied
119    lw_derivatives(icol, nlev+1) = 1.0_jprb
120
121    ! Move up through the atmosphere computing the derivatives at each
122    ! half-level
123    do jlev = nlev,1,-1
124      lw_derivatives_g = lw_derivatives_g * transmittance(:,jlev)
125      lw_derivatives(icol,jlev) = (1.0_jprb - weight) * lw_derivatives(icol,jlev) &
126           &                    + weight * sum(lw_derivatives_g)
127    end do
128
129    if (lhook) call dr_hook('radiation_lw_derivatives:modify_lw_derivatives_ica',1,hook_handle)
130
131  end subroutine modify_lw_derivatives_ica
132
133
134
135  !---------------------------------------------------------------------
136  ! Calculation for solvers involving multiple regions and matrices
137  subroutine calc_lw_derivatives_matrix(ng, nlev, nreg, icol, transmittance, &
138       &                                u_matrix, flux_up_surf, lw_derivatives)
139
140    use parkind1, only           : jprb
141    use yomhook,  only           : lhook, dr_hook
142
143    use radiation_matrix
144
145    implicit none
146
147    ! Inputs
148    integer,    intent(in) :: ng   ! number of spectral intervals
149    integer,    intent(in) :: nlev ! number of levels
150    integer,    intent(in) :: nreg ! number of regions
151    integer,    intent(in) :: icol ! Index of column for output
152    real(jprb), intent(in) :: transmittance(ng,nreg,nreg,nlev)
153    real(jprb), intent(in) :: u_matrix(nreg,nreg,nlev+1) ! Upward overlap matrix
154    real(jprb), intent(in) :: flux_up_surf(ng) ! Upwelling surface spectral flux (W m-2)
155   
156    ! Output
157    real(jprb), intent(out) :: lw_derivatives(:,:) ! dimensioned (ncol,nlev+1)
158
159    ! Rate of change of spectral flux at a given height with respect
160    ! to the surface value
161    real(jprb) :: lw_derivatives_g_reg(ng,nreg)
162
163    integer    :: jlev
164
165    real(jprb) :: hook_handle
166
167    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_matrix',0,hook_handle)
168
169    ! Initialize the derivatives at the surface; the surface is
170    ! treated as a single clear-sky layer so we only need to put
171    ! values in region 1.
172    lw_derivatives_g_reg = 0.0_jprb
173    lw_derivatives_g_reg(:,1) = flux_up_surf / sum(flux_up_surf)
174    lw_derivatives(icol, nlev+1) = 1.0_jprb
175
176    ! Move up through the atmosphere computing the derivatives at each
177    ! half-level
178    do jlev = nlev,1,-1
179      ! Compute effect of overlap at half-level jlev+1, yielding
180      ! derivatives just above that half-level
181      lw_derivatives_g_reg = singlemat_x_vec(ng,ng,nreg,u_matrix(:,:,jlev+1),lw_derivatives_g_reg)
182
183      ! Compute effect of transmittance of layer jlev, yielding
184      ! derivatives just below the half-level above (jlev)
185      lw_derivatives_g_reg = mat_x_vec(ng,ng,nreg,transmittance(:,:,:,jlev),lw_derivatives_g_reg)
186
187      lw_derivatives(icol, jlev) = sum(lw_derivatives_g_reg)
188    end do
189
190    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_matrix',1,hook_handle)
191
192  end subroutine calc_lw_derivatives_matrix
193
194
195  !---------------------------------------------------------------------
196  ! Calculation for solvers involving multiple regions but no 3D
197  ! effects: the difference from calc_lw_derivatives_matrix is that transmittance
198  ! has one less dimensions
199  subroutine calc_lw_derivatives_region(ng, nlev, nreg, icol, transmittance, &
200       &                                u_matrix, flux_up_surf, lw_derivatives)
201
202    use parkind1, only           : jprb
203    use yomhook,  only           : lhook, dr_hook
204
205    use radiation_matrix
206
207    implicit none
208
209    ! Inputs
210    integer,    intent(in) :: ng   ! number of spectral intervals
211    integer,    intent(in) :: nlev ! number of levels
212    integer,    intent(in) :: nreg ! number of regions
213    integer,    intent(in) :: icol ! Index of column for output
214    real(jprb), intent(in) :: transmittance(ng,nreg,nlev)
215    real(jprb), intent(in) :: u_matrix(nreg,nreg,nlev+1) ! Upward overlap matrix
216    real(jprb), intent(in) :: flux_up_surf(ng) ! Upwelling surface spectral flux (W m-2)
217   
218    ! Output
219    real(jprb), intent(out) :: lw_derivatives(:,:) ! dimensioned (ncol,nlev+1)
220
221    ! Rate of change of spectral flux at a given height with respect
222    ! to the surface value
223    real(jprb) :: lw_derivatives_g_reg(ng,nreg)
224
225    integer    :: jlev
226
227    real(jprb) :: hook_handle
228
229    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_region',0,hook_handle)
230
231    ! Initialize the derivatives at the surface; the surface is
232    ! treated as a single clear-sky layer so we only need to put
233    ! values in region 1.
234    lw_derivatives_g_reg = 0.0_jprb
235    lw_derivatives_g_reg(:,1) = flux_up_surf / sum(flux_up_surf)
236    lw_derivatives(icol, nlev+1) = 1.0_jprb
237
238    ! Move up through the atmosphere computing the derivatives at each
239    ! half-level
240    do jlev = nlev,1,-1
241      ! Compute effect of overlap at half-level jlev+1, yielding
242      ! derivatives just above that half-level
243      lw_derivatives_g_reg = singlemat_x_vec(ng,ng,nreg,u_matrix(:,:,jlev+1),lw_derivatives_g_reg)
244
245      ! Compute effect of transmittance of layer jlev, yielding
246      ! derivatives just below the half-level above (jlev)
247      lw_derivatives_g_reg = transmittance(:,:,jlev) * lw_derivatives_g_reg
248
249      lw_derivatives(icol, jlev) = sum(lw_derivatives_g_reg)
250    end do
251
252    if (lhook) call dr_hook('radiation_lw_derivatives:calc_lw_derivatives_region',1,hook_handle)
253
254  end subroutine calc_lw_derivatives_region
255
256
257end module radiation_lw_derivatives
Note: See TracBrowser for help on using the repository browser.