source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_general_cloud_optics_data.F90 @ 5444

Last change on this file since 5444 was 4677, checked in by idelkadi, 16 months ago

Implementation in the LMDZ code of the double call of the ECRAD radiative transfer code to estimate the 3D radiative effect of clouds.

  • This double call of Ecrad is controlled by the ok_3Deffect logic key.
  • If this key is enabled, 2 files of parameter configuration "namelists" for ECRAD are required at runtime: namelist_ecrad and namelist_ecrad_s2.
  • If this key is deactivated, the configuration and initialization part (reading namelist and netcdf files) is performed only once during simulation (1st call to ECRAD). Otherwise, configuration and initialization are performed each time Ecrad is called.
File size: 13.6 KB
Line 
1! radiation_general_cloud_optics_data.F90 - Type to store generalized cloud optical properties
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! License: see the COPYING file for details
15!
16
17module radiation_general_cloud_optics_data
18
19  use parkind1, only : jprb
20
21  implicit none
22
23  public
24
25  !---------------------------------------------------------------------
26  ! This type holds the configuration information to compute optical
27  ! properties for a particular type of cloud or hydrometeor in one of
28  ! the shortwave or longwave
29  type general_cloud_optics_type
30    ! Band-specific (or g-point-specific) values as a look-up table
31    ! versus effective radius dimensioned (nband,n_effective_radius)
32   
33    ! Extinction coefficient per unit mass (m2 kg-1)
34    real(jprb), allocatable, dimension(:,:) :: &
35         &  mass_ext
36   
37    ! Single-scattering albedo and asymmetry factor (dimensionless)
38    real(jprb), allocatable, dimension(:,:) :: &
39         &  ssa, asymmetry
40
41    ! Number of effective radius coefficients, start value and
42    ! interval in look-up table
43    integer    :: n_effective_radius = 0
44    real(jprb) :: effective_radius_0, d_effective_radius
45
46    ! Name of cloud/precip type (e.g. "liquid", "ice", "rain", "snow")
47    ! and the name of the optics scheme.  These two are used to
48    ! generate the name of the data file from which the coefficients
49    ! are read.
50    character(len=511) :: type_name, scheme_name
51   
52    ! Do we use bands or g-points?
53    logical :: use_bands = .false.
54
55   contains
56     procedure :: setup => setup_general_cloud_optics
57     procedure :: add_optical_properties
58
59  end type general_cloud_optics_type
60
61contains
62
63  ! Provides elemental function "delta_eddington"
64#include "radiation_delta_eddington.h"
65
66  !---------------------------------------------------------------------
67  ! Setup cloud optics coefficients by reading them from a file
68  subroutine setup_general_cloud_optics(this, file_name, specdef, &
69       &                                use_bands, use_thick_averaging, &
70       &                                weighting_temperature, &
71       &                                iverbose)
72
73    use yomhook,                       only : lhook, dr_hook
74    use easy_netcdf,                   only : netcdf_file
75    use radiation_spectral_definition, only : spectral_definition_type
76    use radiation_io,                  only : nulout, nulerr, radiation_abort
77
78    class(general_cloud_optics_type), intent(inout)    :: this
79    character(len=*), intent(in)               :: file_name
80    type(spectral_definition_type), intent(in) :: specdef
81    logical, intent(in), optional              :: use_bands, use_thick_averaging
82    real(jprb), intent(in), optional           :: weighting_temperature ! K
83    integer, intent(in), optional              :: iverbose
84   
85    ! Spectral properties read from file, dimensioned (wavenumber,
86    ! n_effective_radius)
87    real(jprb), dimension(:,:), allocatable :: mass_ext, & ! m2 kg-1
88         &                                     ssa, asymmetry
89
90    ! Reflectance of an infinitely thick cloud, needed for thick
91    ! averaging
92    real(jprb), dimension(:,:), allocatable :: ref_inf
93
94    ! Coordinate variables from file
95    real(jprb), dimension(:), allocatable :: wavenumber       ! cm-1
96    real(jprb), dimension(:), allocatable :: effective_radius ! m
97
98    ! Matrix mapping optical properties in the file to values per
99    ! g-point or band, such that in the thin-averaging case,
100    ! this%mass_ext=matmul(mapping,file%mass_ext), so mapping is
101    ! dimensioned (ngpoint,nwav)
102    real(jprb), dimension(:,:), allocatable :: mapping
103
104    ! The NetCDF file containing the coefficients
105    type(netcdf_file)  :: file
106
107    real(jprb) :: diff_spread
108    integer    :: iverb
109    integer    :: nre  ! Number of effective radii
110    integer    :: nwav ! Number of wavenumbers describing cloud
111
112    logical    :: use_bands_local, use_thick_averaging_local
113
114    real(jprb) :: hook_handle
115
116    if (lhook) call dr_hook('radiation_general_cloud_optics_data:setup',0,hook_handle)
117
118    ! Set local values of optional inputs
119    if (present(iverbose)) then
120      iverb = iverbose
121    else
122      iverb = 2
123    end if
124
125    if (present(use_bands)) then
126      use_bands_local = use_bands
127    else
128      use_bands_local = .false.
129    end if
130
131    if (present(use_thick_averaging)) then
132      use_thick_averaging_local = use_thick_averaging
133    else
134      use_thick_averaging_local = .false.
135    end if
136
137    ! Open the scattering file and configure the way it is read
138    call file%open(trim(file_name), iverbose=iverb)
139    !call file%transpose_matrices()
140
141    ! Read coordinate variables
142    call file%get('wavenumber', wavenumber)
143    call file%get('effective_radius', effective_radius)
144
145    ! Read the band-specific coefficients
146    call file%get('mass_extinction_coefficient', mass_ext)
147    call file%get('single_scattering_albedo', ssa)
148    call file%get('asymmetry_factor', asymmetry)
149
150    ! Close scattering file
151    call file%close()
152
153    ! Check effective radius is evenly spaced
154    nre = size(effective_radius)
155    ! Fractional range of differences, should be near zero for evenly
156    ! spaced data
157    diff_spread = (maxval(effective_radius(2:nre)-effective_radius(1:nre-1))  &
158         &        -minval(effective_radius(2:nre)-effective_radius(1:nre-1))) &
159         &      /  minval(abs(effective_radius(2:nre)-effective_radius(1:nre-1)))
160    if (diff_spread > 0.01_jprb) then
161      write(nulerr, '(a,a,a)') '*** Error: effective_radius in ', &
162           &  trim(file_name), ', is not evenly spaced to 1%'
163      call radiation_abort('Radiation configuration error')
164    end if
165
166    ! Set up effective radius coordinate variable
167    this%n_effective_radius = nre
168    this%effective_radius_0 = effective_radius(1)
169    this%d_effective_radius = effective_radius(2) - effective_radius(1)
170
171    ! Set up weighting
172    if (.not. present(weighting_temperature)) then
173      write(nulerr, '(a)') '*** Error: weighting_temperature not provided'
174      call radiation_abort('Radiation configuration error')
175    end if
176
177    nwav = size(wavenumber)
178
179    ! Define the mapping matrix
180    call specdef%calc_mapping(weighting_temperature, &
181         &                    wavenumber, mapping, use_bands=use_bands)
182
183    ! Thick averaging should be performed on delta-Eddington scaled
184    ! quantities (it makes no difference to thin averaging)
185    call delta_eddington(mass_ext, ssa, asymmetry)
186
187
188    ! Thin averaging
189   ! AI juillet 2023
190    allocate(this%mass_ext(nre,nwav))
191    this%mass_ext  = matmul(mapping, mass_ext)
192    this%ssa       = matmul(mapping, mass_ext*ssa) / this%mass_ext
193    this%asymmetry = matmul(mapping, mass_ext*ssa*asymmetry) / (this%mass_ext*this%ssa)
194   
195    if (use_thick_averaging_local) then
196      ! Thick averaging as described by Edwards and Slingo (1996),
197      ! modifying only the single-scattering albedo
198      allocate(ref_inf(nwav, nre))
199
200      ! Eqs. 18 and 17 of Edwards & Slingo (1996)
201      ref_inf = sqrt((1.0_jprb - ssa) / (1.0_jprb - ssa*asymmetry))
202      ref_inf = (1.0_jprb - ref_inf) / (1.0_jprb + ref_inf)
203      ! Here the left-hand side is actually the averaged ref_inf
204      this%ssa = matmul(mapping, ref_inf)
205      ! Eq. 19 of Edwards and Slingo (1996)
206      this%ssa = 4.0_jprb * this%ssa / ((1.0_jprb + this%ssa)**2 &
207           &  - this%asymmetry * (1.0_jprb - this%ssa)**2)
208
209      deallocate(ref_inf)
210    end if
211
212    deallocate(mapping)
213
214    ! Revert back to unscaled quantities
215    call revert_delta_eddington(this%mass_ext, this%ssa, this%asymmetry)
216
217    if (iverb >= 2) then
218      write(nulout,'(a,a)') '  File: ', trim(file_name)
219      write(nulout,'(a,f7.1,a)') '  Weighting temperature: ', weighting_temperature, ' K'
220      if (use_thick_averaging_local) then
221        write(nulout,'(a)') '  SSA averaging: optically thick limit'
222      else
223        write(nulout,'(a)') '  SSA averaging: optically thin limit'
224      end if
225      if (use_bands_local) then
226        write(nulout,'(a,i0,a)') '  Spectral discretization: ', specdef%nband, ' bands'
227      else
228        write(nulout,'(a,i0,a)') '  Spectral discretization: ', specdef%ng, ' g-points'
229      end if
230      write(nulout,'(a,i0,a,f6.1,a,f6.1,a)') '  Effective radius look-up: ', nre, ' points in range ', &
231           &  effective_radius(1)*1.0e6_jprb, '-', effective_radius(nre)*1.0e6_jprb, ' um'
232      write(nulout,'(a,i0,a,i0,a)') '  Wavenumber range: ', int(specdef%min_wavenumber()), '-', &
233           &  int(specdef%max_wavenumber()), ' cm-1'
234    end if
235
236    if (lhook) call dr_hook('radiation_general_cloud_optics_data:setup',1,hook_handle)
237
238  end subroutine setup_general_cloud_optics
239
240
241  !---------------------------------------------------------------------
242  ! Add the optical properties of a particular cloud type to the
243  ! accumulated optical properties of all cloud types
244  subroutine add_optical_properties(this, ng, nlev, ncol, &
245       &                            cloud_fraction, &
246       &                            water_path, effective_radius, &
247       &                            od, scat_od, scat_asymmetry)
248
249    use yomhook, only : lhook, dr_hook
250
251    class(general_cloud_optics_type), intent(in) :: this
252
253    ! Number of g points, levels and columns
254    integer, intent(in) :: ng, nlev, ncol
255
256    ! Properties of present cloud type, dimensioned (ncol,nlev)
257    real(jprb), intent(in) :: cloud_fraction(:,:)
258    real(jprb), intent(in) :: water_path(:,:)       ! kg m-2
259    real(jprb), intent(in) :: effective_radius(:,:) ! m
260
261    ! Optical properties which are additive per cloud type,
262    ! dimensioned (ng,nlev,ncol)
263    real(jprb), intent(inout), dimension(ng,nlev,ncol) &
264         &  :: od             ! Optical depth of layer
265    real(jprb), intent(inout), dimension(ng,nlev,ncol), optional &
266         &  :: scat_od, &     ! Scattering optical depth of layer
267         &     scat_asymmetry ! Scattering optical depth x asymmetry factor
268
269    real(jprb) :: od_local(ng)
270
271    real(jprb) :: re_index, weight1, weight2
272    integer :: ire
273
274    integer :: jcol, jlev
275
276    real(jprb) :: hook_handle
277
278    if (lhook) call dr_hook('radiation_general_cloud_optics_data:add_optical_properties',0,hook_handle)
279
280    if (present(scat_od)) then
281      do jcol = 1,ncol
282        do jlev = 1,nlev
283          if (cloud_fraction(jcol, jlev) > 0.0_jprb) then
284            re_index = max(1.0_jprb, min(1.0_jprb + (effective_radius(jcol,jlev)-this%effective_radius_0) &
285                 &              / this%d_effective_radius, this%n_effective_radius-0.0001_jprb))
286            ire = int(re_index)
287            weight2 = re_index - ire
288            weight1 = 1.0_jprb - weight2
289            od_local = water_path(jcol, jlev) * (weight1*this%mass_ext(:,ire) &
290                 &                              +weight2*this%mass_ext(:,ire+1))
291            od(:,jlev,jcol) = od(:,jlev,jcol) + od_local
292            od_local = od_local * (weight1*this%ssa(:,ire) &
293                 &                +weight2*this%ssa(:,ire+1))
294            scat_od(:,jlev,jcol) = scat_od(:,jlev,jcol) + od_local
295            scat_asymmetry(:,jlev,jcol) = scat_asymmetry(:,jlev,jcol) &
296                 & + od_local * (weight1*this%asymmetry(:,ire) &
297                 &              +weight2*this%asymmetry(:,ire+1))
298          end if
299        end do
300      end do
301    else
302      ! No scattering: return the absorption optical depth
303      do jcol = 1,ncol
304        do jlev = 1,nlev
305          if (water_path(jcol, jlev) > 0.0_jprb) then
306            re_index = max(1.0, min(1.0_jprb + (effective_radius(jcol,jlev)-this%effective_radius_0) &
307                 &              / this%d_effective_radius, this%n_effective_radius-0.0001_jprb))
308            ire = int(re_index)
309            weight2 = re_index - ire
310            weight1 = 1.0_jprb - weight2
311            od(:,jlev,jcol) = od(:,jlev,jcol) &
312                 &  + water_path(jcol, jlev) * (weight1*this%mass_ext(:,ire) &
313                 &                             +weight2*this%mass_ext(:,ire+1)) &
314                 &  * (1.0_jprb - (weight1*this%ssa(:,ire)+weight2*this%ssa(:,ire+1)))
315          end if
316        end do
317      end do
318    end if
319
320    if (lhook) call dr_hook('radiation_general_cloud_optics_data:add_optical_properties',1,hook_handle)
321
322  end subroutine add_optical_properties
323
324
325  !---------------------------------------------------------------------
326  ! Return the Planck function (in W m-2 (cm-1)-1) for a given
327  ! wavenumber (cm-1) and temperature (K), ensuring double precision
328  ! for internal calculation
329  elemental function calc_planck_function_wavenumber(wavenumber, temperature)
330
331    use parkind1,            only : jprb, jprd
332    use radiation_constants, only : SpeedOfLight, BoltzmannConstant, PlanckConstant
333
334    real(jprb), intent(in) :: wavenumber  ! cm-1
335    real(jprb), intent(in) :: temperature ! K
336    real(jprb) :: calc_planck_function_wavenumber
337
338    real(jprd) :: freq ! Hz
339    real(jprd) :: planck_fn_freq ! W m-2 Hz-1
340
341    freq = 100.0_jprd * real(SpeedOfLight,jprd) * real(wavenumber,jprd)
342    planck_fn_freq = 2.0_jprd * real(PlanckConstant,jprd) * freq**3 &
343         &  / (real(SpeedOfLight,jprd)**2 * (exp(real(PlanckConstant,jprd)*freq &
344         &     /(real(BoltzmannConstant,jprd)*real(temperature,jprd))) - 1.0_jprd))
345    calc_planck_function_wavenumber = real(planck_fn_freq * 100.0_jprd * real(SpeedOfLight,jprd), jprb)
346
347  end function calc_planck_function_wavenumber
348
349end module radiation_general_cloud_optics_data
Note: See TracBrowser for help on using the repository browser.