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

Last change on this file since 4582 was 4489, checked in by lguez, 21 months ago

Merge LMDZ_ECRad branch back into trunk!

File size: 11.2 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      allocate(config%cloud_optics_sw(config%n_cloud_types))
78    end if
79
80    if (config%do_lw) then
81      allocate(config%cloud_optics_lw(config%n_cloud_types))
82    end if
83
84    ! Load cloud optics data
85    do jtype = 1,config%n_cloud_types
86      if (config%cloud_type_name(jtype)(1:1) == '/') then
87        file_name = trim(config%cloud_type_name(jtype))
88      else
89        strlen = len_trim(config%cloud_type_name(jtype))
90        if (config%cloud_type_name(jtype)(strlen-2:strlen) == ".nc") then
91          file_name = trim(config%directory_name) &
92               &  // '/' // trim(config%cloud_type_name(jtype))
93        else
94          file_name = trim(config%directory_name) &
95               &  // '/' // trim(config%cloud_type_name(jtype)) &
96               &  // '_scattering.nc'
97        end if
98      end if
99
100      if (config%do_sw) then
101        if (config%iverbosesetup >= 2) then
102          write(nulout,'(a,i0,a)') 'Shortwave cloud type ', jtype, ':'
103        end if
104        call config%cloud_optics_sw(jtype)%setup(file_name, &
105             &  config%gas_optics_sw%spectral_def, &
106             &  use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point), &
107             &  use_thick_averaging=config%use_thick_cloud_spectral_averaging(jtype), &
108             &  weighting_temperature=SolarReferenceTemperature, &
109             &  iverbose=config%iverbosesetup)
110      end if
111
112      if (config%do_lw) then
113        if (config%iverbosesetup >= 2) then
114          write(nulout,'(a,i0,a)') 'Longwave cloud type ', jtype, ':'
115        end if
116        call config%cloud_optics_lw(jtype)%setup(file_name, &
117             &  config%gas_optics_lw%spectral_def, &
118             &  use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point), &
119             &  use_thick_averaging=config%use_thick_cloud_spectral_averaging(jtype), &
120             &  weighting_temperature=TerrestrialReferenceTemperature, &
121             &  iverbose=config%iverbosesetup)
122      end if
123
124    end do
125
126    if (lhook) call dr_hook('radiation_general_cloud_optics:setup_general_cloud_optics',1,hook_handle)
127
128  end subroutine setup_general_cloud_optics
129
130  !---------------------------------------------------------------------
131  ! Compute cloud optical properties
132  subroutine general_cloud_optics(nlev,istartcol,iendcol, &
133       &  config, thermodynamics, cloud, &
134       &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud, &
135       &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
136
137    use parkind1, only           : jprb
138    use yomhook,  only           : lhook, dr_hook
139
140    use radiation_io,     only : nulout
141    use radiation_config, only : config_type
142    use radiation_thermodynamics, only    : thermodynamics_type
143    use radiation_cloud, only             : cloud_type
144    use radiation_constants, only         : AccelDueToGravity
145    !use radiation_general_cloud_optics_data, only : general_cloud_optics_type
146
147    integer, intent(in) :: nlev               ! number of model levels
148    integer, intent(in) :: istartcol, iendcol ! range of columns to process
149    type(config_type), intent(in), target :: config
150    type(thermodynamics_type),intent(in)  :: thermodynamics
151    type(cloud_type),   intent(in)        :: cloud
152
153    ! Layer optical depth, single scattering albedo and g factor of
154    ! clouds in each longwave band, where the latter two
155    ! variables are only defined if cloud longwave scattering is
156    ! enabled (otherwise both are treated as zero).
157    real(jprb), dimension(config%n_bands_lw,nlev,istartcol:iendcol), intent(out) :: &
158         &   od_lw_cloud
159    real(jprb), dimension(config%n_bands_lw_if_scattering,nlev,istartcol:iendcol), &
160         &   intent(out) :: ssa_lw_cloud, g_lw_cloud
161
162    ! Layer optical depth, single scattering albedo and g factor of
163    ! clouds in each shortwave band
164    real(jprb), dimension(config%n_bands_sw,nlev,istartcol:iendcol), intent(out) :: &
165         &   od_sw_cloud, ssa_sw_cloud, g_sw_cloud
166
167    ! In-cloud water path of one cloud type (kg m-2)
168    real(jprb), dimension(istartcol:iendcol,nlev) :: water_path
169
170    integer :: jtype, jcol, jlev
171
172    real(jprb) :: hook_handle
173
174    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',0,hook_handle)
175
176    if (config%iverbose >= 2) then
177      write(nulout,'(a)') 'Computing cloud absorption/scattering properties'
178    end if
179
180    ! Array-wise assignment
181    od_lw_cloud  = 0.0_jprb
182    od_sw_cloud  = 0.0_jprb
183    ssa_sw_cloud = 0.0_jprb
184    g_sw_cloud   = 0.0_jprb
185    if (config%do_lw_cloud_scattering) then
186      ssa_lw_cloud = 0.0_jprb
187      g_lw_cloud   = 0.0_jprb
188    end if
189
190    ! Loop over cloud types
191    do jtype = 1,config%n_cloud_types
192      ! Compute in-cloud water path
193      if (config%is_homogeneous) then
194        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
195             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
196             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
197             &  * (1.0_jprb / AccelDueToGravity)
198      else
199        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
200             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
201             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
202             &  * (1.0_jprb / (AccelDueToGravity &
203             &                 * max(config%cloud_fraction_threshold, &
204             &                       cloud%fraction(istartcol:iendcol,:))))
205      end if
206
207      ! Add optical properties to the cumulative total for the
208      ! longwave and shortwave
209      if (config%do_lw) then
210        ! For the moment, we use ssa_lw_cloud and g_lw_cloud as
211        ! containers for scattering optical depth and scattering
212        ! coefficient x asymmetry factor, then scale after
213        if (config%do_lw_cloud_scattering) then
214          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
215               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
216               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
217               &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud)
218        else
219          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
220               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
221               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), od_lw_cloud)
222        end if
223      end if
224     
225      if (config%do_sw) then
226        ! For the moment, we use ssa_sw_cloud and g_sw_cloud as
227        ! containers for scattering optical depth and scattering
228        ! coefficient x asymmetry factor, then scale after
229        call config%cloud_optics_sw(jtype)%add_optical_properties(config%n_bands_sw, nlev, &
230             &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
231             &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
232             &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
233      end if
234    end do
235
236    ! Scale the combined longwave optical properties
237    if (config%do_lw_cloud_scattering) then
238      do jcol = istartcol, iendcol
239        do jlev = 1,nlev
240          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
241            ! Note that original cloud optics does not do
242            ! delta-Eddington scaling for liquid clouds in longwave
243            call delta_eddington_extensive(od_lw_cloud(:,jlev,jcol), &
244                 &  ssa_lw_cloud(:,jlev,jcol), g_lw_cloud(:,jlev,jcol))
245           
246            ! Scale to get asymmetry factor and single scattering albedo
247            g_lw_cloud(:,jlev,jcol) = g_lw_cloud(:,jlev,jcol) &
248                 &  / max(ssa_lw_cloud(:,jlev,jcol), 1.0e-15_jprb)
249            ssa_lw_cloud(:,jlev,jcol) = ssa_lw_cloud(:,jlev,jcol) &
250                 &  / max(od_lw_cloud(:,jlev,jcol),  1.0e-15_jprb)
251          end if
252        end do
253      end do
254    end if
255   
256    ! Scale the combined shortwave optical properties
257    if (config%do_sw) then
258      if (.not. config%do_sw_delta_scaling_with_gases) then
259        do jcol = istartcol, iendcol
260          do jlev = 1,nlev
261            if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
262              call delta_eddington_extensive(od_sw_cloud(:,jlev,jcol), &
263                   &  ssa_sw_cloud(:,jlev,jcol), g_sw_cloud(:,jlev,jcol))
264            end if
265          end do
266        end do
267      end if
268
269      do jcol = istartcol, iendcol
270        do jlev = 1,nlev
271          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
272            ! Scale to get asymmetry factor and single scattering albedo
273            g_sw_cloud(:,jlev,jcol) = g_sw_cloud(:,jlev,jcol) &
274                 &  / max(ssa_sw_cloud(:,jlev,jcol), 1.0e-15_jprb)
275            ssa_sw_cloud(:,jlev,jcol) = ssa_sw_cloud(:,jlev,jcol) &
276                 &  / max(od_sw_cloud(:,jlev,jcol),  1.0e-15_jprb)
277          end if
278        end do
279      end do
280    end if
281
282    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',1,hook_handle)
283
284  end subroutine general_cloud_optics
285
286end module radiation_general_cloud_optics
Note: See TracBrowser for help on using the repository browser.