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

Last change on this file since 4773 was 4773, checked in by idelkadi, 5 months ago
  • Update of Ecrad in LMDZ The same organization of the Ecrad offline version is retained in order to facilitate the updating of Ecrad in LMDZ and the comparison between online and offline results. version 1.6.1 of Ecrad (https://github.com/lguez/ecrad.git)
  • Implementation of the double call of Ecrad in LMDZ


File size: 12.8 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, jphook
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(jphook) :: 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        config%cloud_optics_sw(jtype)%type_name = trim(config%cloud_type_name(jtype))
113      end if
114
115      if (config%do_lw) then
116        if (config%iverbosesetup >= 2) then
117          write(nulout,'(a,i0,a)') 'Longwave cloud type ', jtype, ':'
118        end if
119        call config%cloud_optics_lw(jtype)%setup(file_name, &
120             &  config%gas_optics_lw%spectral_def, &
121             &  use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point), &
122             &  use_thick_averaging=config%use_thick_cloud_spectral_averaging(jtype), &
123             &  weighting_temperature=TerrestrialReferenceTemperature, &
124             &  iverbose=config%iverbosesetup)
125        config%cloud_optics_lw(jtype)%type_name = trim(config%cloud_type_name(jtype))
126      end if
127
128    end do
129
130    if (lhook) call dr_hook('radiation_general_cloud_optics:setup_general_cloud_optics',1,hook_handle)
131
132  end subroutine setup_general_cloud_optics
133
134  !---------------------------------------------------------------------
135  ! Compute cloud optical properties
136  subroutine general_cloud_optics(nlev,istartcol,iendcol, &
137       &  config, thermodynamics, cloud, &
138       &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud, &
139       &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
140
141    use parkind1, only           : jprb
142    use yomhook,  only           : lhook, dr_hook, jphook
143
144    use radiation_io,     only : nulout
145    use radiation_config, only : config_type
146    use radiation_thermodynamics, only    : thermodynamics_type
147    use radiation_cloud, only             : cloud_type
148    use radiation_constants, only         : AccelDueToGravity
149    !use radiation_general_cloud_optics_data, only : general_cloud_optics_type
150
151    integer, intent(in) :: nlev               ! number of model levels
152    integer, intent(in) :: istartcol, iendcol ! range of columns to process
153    type(config_type), intent(in), target :: config
154    type(thermodynamics_type),intent(in)  :: thermodynamics
155    type(cloud_type),   intent(in)        :: cloud
156
157    ! Layer optical depth, single scattering albedo and g factor of
158    ! clouds in each longwave band, where the latter two
159    ! variables are only defined if cloud longwave scattering is
160    ! enabled (otherwise both are treated as zero).
161    real(jprb), dimension(config%n_bands_lw,nlev,istartcol:iendcol), intent(out) :: &
162         &   od_lw_cloud
163    real(jprb), dimension(config%n_bands_lw_if_scattering,nlev,istartcol:iendcol), &
164         &   intent(out) :: ssa_lw_cloud, g_lw_cloud
165
166    ! Layer optical depth, single scattering albedo and g factor of
167    ! clouds in each shortwave band
168    real(jprb), dimension(config%n_bands_sw,nlev,istartcol:iendcol), intent(out) :: &
169         &   od_sw_cloud, ssa_sw_cloud, g_sw_cloud
170
171    ! In-cloud water path of one cloud type (kg m-2)
172    real(jprb), dimension(istartcol:iendcol,nlev) :: water_path
173
174    integer :: jtype, jcol, jlev
175
176    real(jphook) :: hook_handle
177
178    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',0,hook_handle)
179
180    if (config%iverbose >= 2) then
181      write(nulout,'(a)') 'Computing cloud absorption/scattering properties'
182    end if
183
184    ! Array-wise assignment
185    od_lw_cloud  = 0.0_jprb
186    od_sw_cloud  = 0.0_jprb
187    ssa_sw_cloud = 0.0_jprb
188    g_sw_cloud   = 0.0_jprb
189    if (config%do_lw_cloud_scattering) then
190      ssa_lw_cloud = 0.0_jprb
191      g_lw_cloud   = 0.0_jprb
192    end if
193
194    ! Loop over cloud types
195    do jtype = 1,config%n_cloud_types
196      ! Compute in-cloud water path
197      if (config%is_homogeneous) then
198        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
199             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
200             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
201             &  * (1.0_jprb / AccelDueToGravity)
202      else
203        water_path = cloud%mixing_ratio(istartcol:iendcol,:,jtype) &
204             &  *  (thermodynamics%pressure_hl(istartcol:iendcol, 2:nlev+1) &
205             &     -thermodynamics%pressure_hl(istartcol:iendcol, 1:nlev)) &
206             &  * (1.0_jprb / (AccelDueToGravity &
207             &                 * max(config%cloud_fraction_threshold, &
208             &                       cloud%fraction(istartcol:iendcol,:))))
209      end if
210
211      ! Add optical properties to the cumulative total for the
212      ! longwave and shortwave
213      if (config%do_lw) then
214        ! For the moment, we use ssa_lw_cloud and g_lw_cloud as
215        ! containers for scattering optical depth and scattering
216        ! coefficient x asymmetry factor, then scale after
217        if (config%do_lw_cloud_scattering) then
218          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
219               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
220               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
221               &  od_lw_cloud, ssa_lw_cloud, g_lw_cloud)
222        else
223          call config%cloud_optics_lw(jtype)%add_optical_properties(config%n_bands_lw, nlev, &
224               &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
225               &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), od_lw_cloud)
226        end if
227      end if
228     
229      if (config%do_sw) then
230        ! For the moment, we use ssa_sw_cloud and g_sw_cloud as
231        ! containers for scattering optical depth and scattering
232        ! coefficient x asymmetry factor, then scale after
233        call config%cloud_optics_sw(jtype)%add_optical_properties(config%n_bands_sw, nlev, &
234             &  iendcol+1-istartcol, cloud%fraction(istartcol:iendcol,:), &
235             &  water_path, cloud%effective_radius(istartcol:iendcol,:,jtype), &
236             &  od_sw_cloud, ssa_sw_cloud, g_sw_cloud)
237      end if
238    end do
239
240    ! Scale the combined longwave optical properties
241    if (config%do_lw_cloud_scattering) then
242      do jcol = istartcol, iendcol
243        do jlev = 1,nlev
244          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
245            ! Note that original cloud optics does not do
246            ! delta-Eddington scaling for liquid clouds in longwave
247            call delta_eddington_extensive(od_lw_cloud(:,jlev,jcol), &
248                 &  ssa_lw_cloud(:,jlev,jcol), g_lw_cloud(:,jlev,jcol))
249           
250            ! Scale to get asymmetry factor and single scattering albedo
251            g_lw_cloud(:,jlev,jcol) = g_lw_cloud(:,jlev,jcol) &
252                 &  / max(ssa_lw_cloud(:,jlev,jcol), 1.0e-15_jprb)
253            ssa_lw_cloud(:,jlev,jcol) = ssa_lw_cloud(:,jlev,jcol) &
254                 &  / max(od_lw_cloud(:,jlev,jcol),  1.0e-15_jprb)
255          end if
256        end do
257      end do
258    end if
259   
260    ! Scale the combined shortwave optical properties
261    if (config%do_sw) then
262      if (.not. config%do_sw_delta_scaling_with_gases) then
263        do jcol = istartcol, iendcol
264          do jlev = 1,nlev
265            if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
266              call delta_eddington_extensive(od_sw_cloud(:,jlev,jcol), &
267                   &  ssa_sw_cloud(:,jlev,jcol), g_sw_cloud(:,jlev,jcol))
268            end if
269          end do
270        end do
271      end if
272
273      do jcol = istartcol, iendcol
274        do jlev = 1,nlev
275          if (cloud%fraction(jcol,jlev) > 0.0_jprb) then
276            ! Scale to get asymmetry factor and single scattering albedo
277            g_sw_cloud(:,jlev,jcol) = g_sw_cloud(:,jlev,jcol) &
278                 &  / max(ssa_sw_cloud(:,jlev,jcol), 1.0e-15_jprb)
279            ssa_sw_cloud(:,jlev,jcol) = ssa_sw_cloud(:,jlev,jcol) &
280                 &  / max(od_sw_cloud(:,jlev,jcol),  1.0e-15_jprb)
281          end if
282        end do
283      end do
284    end if
285
286    if (lhook) call dr_hook('radiation_general_cloud_optics:general_cloud_optics',1,hook_handle)
287
288  end subroutine general_cloud_optics
289
290
291  !---------------------------------------------------------------------
292  ! Save all the cloud optics look-up tables for sw/lw and for each
293  ! hydrometeor type
294  subroutine save_general_cloud_optics(config, file_prefix, iverbose)
295
296    use yomhook,     only : lhook, dr_hook, jphook
297    use easy_netcdf, only : netcdf_file
298    use radiation_config, only : config_type
299   
300    type(config_type),  intent(in) :: config
301    character(len=*),   intent(in) :: file_prefix
302    integer,  optional, intent(in) :: iverbose
303
304    integer :: jtype
305
306    real(jphook) :: hook_handle
307
308    if (lhook) call dr_hook('radiation_general_cloud_optics:save',0,hook_handle)
309
310    do jtype = 1,config%n_cloud_types
311      if (config%do_sw) then
312        associate(co_sw => config%cloud_optics_sw(jtype))
313          call co_sw%save(file_prefix//"_sw_" &
314               &          //trim(co_sw%type_name)//".nc", iverbose)
315        end associate
316      end if
317
318      if (config%do_lw) then
319        associate(co_lw => config%cloud_optics_lw(jtype))
320          call co_lw%save(file_prefix//"_lw_" &
321               &          //trim(co_lw%type_name)//".nc", iverbose)
322        end associate
323      end if
324    end do
325   
326    if (lhook) call dr_hook('radiation_general_cloud_optics:save',1,hook_handle)
327
328  end subroutine save_general_cloud_optics
329
330end module radiation_general_cloud_optics
Note: See TracBrowser for help on using the repository browser.