source: LMDZ6/trunk/libf/phylmd/ecrad/radiation_general_cloud_optics.F90 @ 4677

Last change on this file since 4677 was 4677, checked in by idelkadi, 9 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: 11.4 KB
Line 
1! radiation_general_cloud_optics.F90 - Computing generalized cloud optical properties
2!
3! (C) Copyright 2020- 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
18
19  implicit none
20
21  public
22 
23contains
24
25  ! Provides elemental function "delta_eddington_scat_od"
26#include "radiation_delta_eddington.h"
27
28
29  !---------------------------------------------------------------------
30  ! Load cloud scattering data; this subroutine delegates to one
31  ! in radiation_general_cloud_optics_data.F90
32  subroutine setup_general_cloud_optics(config)
33
34    use parkind1,         only : jprb
35    use yomhook,          only : lhook, dr_hook
36
37    use radiation_io,     only : nulout
38    use radiation_config, only : config_type, NMaxCloudTypes
39    use radiation_spectral_definition, only : SolarReferenceTemperature, &
40         &                                    TerrestrialReferenceTemperature
41
42    type(config_type), intent(inout) :: config
43
44    character(len=511) :: file_name
45
46    integer :: jtype ! loop index
47    integer :: strlen
48
49    real(jprb) :: hook_handle
50
51    if (lhook) call dr_hook('radiation_general_cloud_optics:setup_general_cloud_optics',0,hook_handle)
52
53    ! Count number of cloud types
54    config%n_cloud_types = 0
55    do jtype = 1,NMaxCloudTypes
56      if (len_trim(config%cloud_type_name(jtype)) > 0) then
57        config%n_cloud_types = jtype
58      else
59        exit
60      end if
61    end do
62
63    ! If cloud_type_name has not been provided then assume liquid,ice
64    ! noting that the default spectral averaging (defined in
65    ! radiation_config.F90) is "thick"
66    if (config%n_cloud_types == 0) then
67      config%cloud_type_name(1) = "mie_droplet"
68      config%cloud_type_name(2) = "baum-general-habit-mixture_ice"
69      ! Optionally override spectral averaging method
70      !config%use_thick_cloud_spectral_averaging(1) = .false.
71      !config%use_thick_cloud_spectral_averaging(2) = .false.
72      config%n_cloud_types = 2
73    end if
74
75    ! Allocate structures
76    if (config%do_sw) then
77      if (allocated(config%cloud_optics_sw)) deallocate(config%cloud_optics_sw)     
78      allocate(config%cloud_optics_sw(config%n_cloud_types))
79    end if
80
81    if (config%do_lw) then
82      if (allocated(config%cloud_optics_lw)) deallocate(config%cloud_optics_lw)       
83      allocate(config%cloud_optics_lw(config%n_cloud_types))
84    end if
85
86    ! Load cloud optics data
87    do jtype = 1,config%n_cloud_types
88      if (config%cloud_type_name(jtype)(1:1) == '/') then
89        file_name = trim(config%cloud_type_name(jtype))
90      else
91        strlen = len_trim(config%cloud_type_name(jtype))
92        if (config%cloud_type_name(jtype)(strlen-2:strlen) == ".nc") then
93          file_name = trim(config%directory_name) &
94               &  // '/' // trim(config%cloud_type_name(jtype))
95        else
96          file_name = trim(config%directory_name) &
97               &  // '/' // trim(config%cloud_type_name(jtype)) &
98               &  // '_scattering.nc'
99        end if
100      end if
101
102      if (config%do_sw) then
103        if (config%iverbosesetup >= 2) then
104          write(nulout,'(a,i0,a)') 'Shortwave cloud type ', jtype, ':'
105        end if
106        call config%cloud_optics_sw(jtype)%setup(file_name, &
107             &  config%gas_optics_sw%spectral_def, &
108             &  use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point), &
109             &  use_thick_averaging=config%use_thick_cloud_spectral_averaging(jtype), &
110             &  weighting_temperature=SolarReferenceTemperature, &
111             &  iverbose=config%iverbosesetup)
112      end if
113
114      if (config%do_lw) then
115        if (config%iverbosesetup >= 2) then
116          write(nulout,'(a,i0,a)') 'Longwave cloud type ', jtype, ':'
117        end if
118        call config%cloud_optics_lw(jtype)%setup(file_name, &
119             &  config%gas_optics_lw%spectral_def, &
120             &  use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point), &
121             &  use_thick_averaging=config%use_thick_cloud_spectral_averaging(jtype), &
122             &  weighting_temperature=TerrestrialReferenceTemperature, &
123             &  iverbose=config%iverbosesetup)
124      end if
125
126    end do
127
128    if (lhook) call dr_hook('radiation_general_cloud_optics:setup_general_cloud_optics',1,hook_handle)
129
130  end subroutine setup_general_cloud_optics
131
132  !---------------------------------------------------------------------
133  ! Compute cloud optical properties
134  subroutine general_cloud_optics(nlev,istartcol,iendcol, &
135       &  config, thermodynamics, cloud, &
136       &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud, &
137       &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
138
139    use parkind1, only           : jprb
140    use yomhook,  only           : lhook, dr_hook
141
142    use radiation_io,     only : nulout
143    use radiation_config, only : config_type
144    use radiation_thermodynamics, only    : thermodynamics_type
145    use radiation_cloud, only             : cloud_type
146    use radiation_constants, only         : AccelDueToGravity
147    !use radiation_general_cloud_optics_data, only : general_cloud_optics_type
148
149    integer, intent(in) :: nlev               ! number of model levels
150    integer, intent(in) :: istartcol, iendcol ! range of columns to process
151    type(config_type), intent(in), target :: config
152    type(thermodynamics_type),intent(in)  :: thermodynamics
153    type(cloud_type),   intent(in)        :: cloud
154
155    ! Layer optical depth, single scattering albedo and g factor of
156    ! clouds in each longwave band, where the latter two
157    ! variables are only defined if cloud longwave scattering is
158    ! enabled (otherwise both are treated as zero).
159    real(jprb), dimension(config%n_bands_lw,nlev,istartcol:iendcol), intent(out) :: &
160         &   od_lw_cloud
161    real(jprb), dimension(config%n_bands_lw_if_scattering,nlev,istartcol:iendcol), &
162         &   intent(out) :: ssa_lw_cloud, g_lw_cloud
163
164    ! Layer optical depth, single scattering albedo and g factor of
165    ! clouds in each shortwave band
166    real(jprb), dimension(config%n_bands_sw,nlev,istartcol:iendcol), intent(out) :: &
167         &   od_sw_cloud, ssa_sw_cloud, g_sw_cloud
168
169    ! In-cloud water path of one cloud type (kg m-2)
170    real(jprb), dimension(istartcol:iendcol,nlev) :: water_path
171
172    integer :: jtype, jcol, jlev
173
174    real(jprb) :: hook_handle
175
176    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',0,hook_handle)
177
178    if (config%iverbose >= 2) then
179      write(nulout,'(a)') 'Computing cloud absorption/scattering properties'
180    end if
181
182    ! Array-wise assignment
183    od_lw_cloud  = 0.0_jprb
184    od_sw_cloud  = 0.0_jprb
185    ssa_sw_cloud = 0.0_jprb
186    g_sw_cloud   = 0.0_jprb
187    if (config%do_lw_cloud_scattering) then
188      ssa_lw_cloud = 0.0_jprb
189      g_lw_cloud   = 0.0_jprb
190    end if
191
192    ! Loop over cloud types
193    do jtype = 1,config%n_cloud_types
194      ! Compute in-cloud water path
195      if (config%is_homogeneous) then
196        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
197             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
198             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
199             &  * (1.0_jprb / AccelDueToGravity)
200      else
201        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
202             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
203             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
204             &  * (1.0_jprb / (AccelDueToGravity &
205             &                 * max(config%cloud_fraction_threshold, &
206             &                       cloud%fraction(istartcol:iendcol,:))))
207      end if
208
209      ! Add optical properties to the cumulative total for the
210      ! longwave and shortwave
211      if (config%do_lw) then
212        ! For the moment, we use ssa_lw_cloud and g_lw_cloud as
213        ! containers for scattering optical depth and scattering
214        ! coefficient x asymmetry factor, then scale after
215        if (config%do_lw_cloud_scattering) then
216          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
217               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
218               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
219               &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud)
220        else
221          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
222               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
223               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), od_lw_cloud)
224        end if
225      end if
226     
227      if (config%do_sw) then
228        ! For the moment, we use ssa_sw_cloud and g_sw_cloud as
229        ! containers for scattering optical depth and scattering
230        ! coefficient x asymmetry factor, then scale after
231        call config%cloud_optics_sw(jtype)%add_optical_properties(config%n_bands_sw, nlev, &
232             &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
233             &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
234             &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
235      end if
236    end do
237
238    ! Scale the combined longwave optical properties
239    if (config%do_lw_cloud_scattering) then
240      do jcol = istartcol, iendcol
241        do jlev = 1,nlev
242          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
243            ! Note that original cloud optics does not do
244            ! delta-Eddington scaling for liquid clouds in longwave
245            call delta_eddington_extensive(od_lw_cloud(:,jlev,jcol), &
246                 &  ssa_lw_cloud(:,jlev,jcol), g_lw_cloud(:,jlev,jcol))
247           
248            ! Scale to get asymmetry factor and single scattering albedo
249            g_lw_cloud(:,jlev,jcol) = g_lw_cloud(:,jlev,jcol) &
250                 &  / max(ssa_lw_cloud(:,jlev,jcol), 1.0e-15_jprb)
251            ssa_lw_cloud(:,jlev,jcol) = ssa_lw_cloud(:,jlev,jcol) &
252                 &  / max(od_lw_cloud(:,jlev,jcol),  1.0e-15_jprb)
253          end if
254        end do
255      end do
256    end if
257   
258    ! Scale the combined shortwave optical properties
259    if (config%do_sw) then
260      if (.not. config%do_sw_delta_scaling_with_gases) then
261        do jcol = istartcol, iendcol
262          do jlev = 1,nlev
263            if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
264              call delta_eddington_extensive(od_sw_cloud(:,jlev,jcol), &
265                   &  ssa_sw_cloud(:,jlev,jcol), g_sw_cloud(:,jlev,jcol))
266            end if
267          end do
268        end do
269      end if
270
271      do jcol = istartcol, iendcol
272        do jlev = 1,nlev
273          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
274            ! Scale to get asymmetry factor and single scattering albedo
275            g_sw_cloud(:,jlev,jcol) = g_sw_cloud(:,jlev,jcol) &
276                 &  / max(ssa_sw_cloud(:,jlev,jcol), 1.0e-15_jprb)
277            ssa_sw_cloud(:,jlev,jcol) = ssa_sw_cloud(:,jlev,jcol) &
278                 &  / max(od_sw_cloud(:,jlev,jcol),  1.0e-15_jprb)
279          end if
280        end do
281      end do
282    end if
283
284    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',1,hook_handle)
285
286  end subroutine general_cloud_optics
287
288end module radiation_general_cloud_optics
Note: See TracBrowser for help on using the repository browser.