source: LMDZ6/trunk/libf/phylmd/ecrad.v1.5.1/radiation_aerosol_optics.F90 @ 5213

Last change on this file since 5213 was 4518, checked in by lguez, 19 months ago

Remove setup_general_aerosol_optics_lmdz

Remove procedure setup_general_aerosol_optics_lmdz. This procedure
was added in commit r4444. For the now moment, we do not need a
specific procedure for the case config%use_general_aerosol_optics :
"Read file containing high spectral resolution optical properties and
average to the spectral intervals of the current gas-optics scheme".

File size: 46.2 KB
Line 
1! radiation_aerosol_optics.F90 - Computing aerosol optical properties
2!
3! (C) Copyright 2015- 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! Modifications
16!   2018-04-15  R. Hogan  Add "direct" option
17!   2020-11-14  R. Hogan  Add setup_general_aerosol_optics for ecCKD compatibility
18!   2022-03-27  R. Hogan  Add setup_general_aerosol_optics_legacy to use RRTM aerosol files with ecCKD
19
20module radiation_aerosol_optics
21
22  implicit none
23  public
24
25contains
26
27  ! Provides the elemental function "delta_eddington_extensive"
28#include "radiation_delta_eddington.h"
29
30  !---------------------------------------------------------------------
31  ! Load aerosol scattering data; this subroutine delegates to one
32  ! in radiation_aerosol_optics_data.F90
33  subroutine setup_aerosol_optics(config)
34
35    use parkind1,                      only : jprb
36    use yomhook,                       only : lhook, dr_hook
37    use radiation_config,              only : config_type
38    use radiation_aerosol_optics_data, only : aerosol_optics_type
39    use radiation_io,                  only : nulerr, radiation_abort
40    use setup_aerosol_optics_lmdz_m,   only: setup_aerosol_optics_lmdz
41
42    type(config_type), intent(inout) :: config
43
44    real(jprb) :: hook_handle
45
46    if (lhook) call dr_hook('radiation_aerosol_optics:setup_aerosol_optics',0,hook_handle)
47
48    if (config%n_aerosol_types > 0) then
49      ! Load data from file and prepare to map config%n_aerosol_types
50      ! aerosol types
51      if (config%use_general_aerosol_optics) then
52        ! Read file containing high spectral resolution optical
53        ! properties and average to the spectral intervals of the
54        ! current gas-optics scheme
55        call setup_general_aerosol_optics(config)
56      else
57        ! Read file containing optical properties already in the bands
58        ! of the gas-optics scheme
59!        call config%aerosol_optics%setup(trim(config%aerosol_optics_file_name), &
60!             &                           iverbose=config%iverbosesetup)
61        call setup_aerosol_optics_lmdz(config%aerosol_optics, &
62                                       trim(config%aerosol_optics_file_name))
63      end if
64
65      call config%aerosol_optics%initialize_types(config%n_aerosol_types)
66
67      ! Check agreement in number of bands
68      if (config%n_bands_lw /= config%aerosol_optics%n_bands_lw) then
69        write(nulerr,'(a,i0,a,i0,a)') '*** Error: number of longwave bands (', &
70             &  config%n_bands_lw, ') does not match aerosol optics look-up table (', &
71             &  config%aerosol_optics%n_bands_lw, ')'
72        call radiation_abort()
73      end if
74      if (config%n_bands_sw /= config%aerosol_optics%n_bands_sw) then
75        write(nulerr,'(a)') '*** Error: number of shortwave bands does not match aerosol optics look-up table'
76        call radiation_abort()
77      end if
78
79      ! Map aerosol types to those loaded from the data file
80      call config%aerosol_optics%set_types(config%i_aerosol_type_map(1:config%n_aerosol_types))
81    end if
82
83    if (config%iverbosesetup >= 1) then
84      call config%aerosol_optics%print_description(config%i_aerosol_type_map(1:config%n_aerosol_types))
85    end if
86
87    if (lhook) call dr_hook('radiation_aerosol_optics:setup_aerosol_optics',1,hook_handle)
88
89  end subroutine setup_aerosol_optics
90
91
92  !---------------------------------------------------------------------
93  ! Read file containing high spectral resolution optical properties
94  ! and average to the spectral intervals of the current gas-optics
95  ! scheme
96  subroutine setup_general_aerosol_optics(config)
97
98    use parkind1,                      only : jprb
99    use yomhook,                       only : lhook, dr_hook
100    use easy_netcdf,                   only : netcdf_file
101    use radiation_config,              only : config_type
102    use radiation_aerosol_optics_data, only : aerosol_optics_type
103    use radiation_spectral_definition, only : SolarReferenceTemperature, &
104         &                                    TerrestrialReferenceTemperature
105    use radiation_io,                  only : nulout
106
107    type(config_type), intent(inout), target :: config
108
109    ! The NetCDF file containing the aerosol optics data
110    type(netcdf_file)  :: file
111
112    ! Wavenumber points in NetCDF file
113    real(jprb), allocatable :: wavenumber(:) ! cm-1
114
115    ! Hydrophilic aerosol properties
116    real(jprb), allocatable :: mass_ext_philic(:,:,:)    ! Mass-ext coefficient (m2 kg-1)
117    real(jprb), allocatable :: ssa_philic(:,:,:)         ! Single-scattering albedo
118    real(jprb), allocatable :: g_philic(:,:,:)           ! Asymmetry factor
119    real(jprb), allocatable :: lidar_ratio_philic(:,:,:) ! Lidar ratio (sr)
120
121    ! Hydrophobic aerosol properties
122    real(jprb), allocatable :: mass_ext_phobic(:,:)      ! Mass-ext coefficient (m2 kg-1)
123    real(jprb), allocatable :: ssa_phobic(:,:)           ! Single-scattering albedo
124    real(jprb), allocatable :: g_phobic(:,:)             ! Asymmetry factor
125    real(jprb), allocatable :: lidar_ratio_phobic(:,:)   ! Lidar ratio (sr)
126
127    ! Mapping matrix between optical properties at the wavenumbers in
128    ! the file, and spectral intervals used by the gas-optics scheme
129    real(jprb), allocatable :: mapping(:,:)
130
131    ! Pointer to the aerosol optics coefficients for brevity of access
132    type(aerosol_optics_type), pointer :: ao
133
134    ! Target monochromatic wavenumber for interpolation (cm-1)
135    real(jprb) :: wavenumber_target
136
137    ! Number of spectral points describing aerosol properties in the
138    ! shortwave and longwave
139    integer    :: nspecsw, nspeclw
140
141    ! Number of monochromatic wavelengths required
142    integer    :: nmono
143
144    integer    :: n_type_philic, n_type_phobic, nrh, nwn
145    integer    :: jtype, jwl, iwn
146
147    ! Weight of first point in interpolation
148    real(jprb) :: weight1
149
150    real(jprb) :: hook_handle
151
152    if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',0,hook_handle)
153
154    ao => config%aerosol_optics
155
156    call file%open(trim(config%aerosol_optics_file_name), iverbose=config%iverbosesetup)
157
158    if (.not. file%exists('wavenumber')) then
159      ! Assume we have an old-style aerosol optics file with optical
160      ! properties provided per pre-defined band
161      call file%close()
162      if (config%iverbosesetup >= 2) then
163        write(nulout,'(a)') 'Legacy aerosol optics file: mapping between bands'
164      end if
165      call setup_general_aerosol_optics_legacy(config, trim(config%aerosol_optics_file_name))
166      return
167    end if
168
169    if (file%exists('mass_ext_hydrophilic')) then
170      ao%use_hydrophilic = .true.
171    else
172      ao%use_hydrophilic = .false.
173    end if
174 
175    call file%get('wavenumber', wavenumber)
176    nwn = size(wavenumber)
177
178    ! Read the raw scattering data
179    call file%get('mass_ext_hydrophobic',    mass_ext_phobic)
180    call file%get('ssa_hydrophobic',         ssa_phobic)
181    call file%get('asymmetry_hydrophobic',   g_phobic)
182    call file%get('lidar_ratio_hydrophobic', lidar_ratio_phobic)
183
184    call file%get_global_attribute('description_hydrophobic', &
185         &                         ao%description_phobic_str)
186
187    if (ao%use_hydrophilic) then
188      call file%get('mass_ext_hydrophilic',    mass_ext_philic)
189      call file%get('ssa_hydrophilic',         ssa_philic)
190      call file%get('asymmetry_hydrophilic',   g_philic)
191      call file%get('lidar_ratio_hydrophilic', lidar_ratio_philic)
192
193      call file%get('relative_humidity1',      ao%rh_lower)
194
195      call file%get_global_attribute('description_hydrophilic', &
196           &                         ao%description_philic_str)
197    end if
198
199    ! Close aerosol scattering file
200    call file%close()
201
202    n_type_phobic = size(mass_ext_phobic, 2)
203    if (ao%use_hydrophilic) then
204      n_type_philic = size(mass_ext_philic, 3)
205      nrh = size(ao%rh_lower)
206    else
207      n_type_philic = 0
208      nrh = 0
209    end if
210
211    if (config%do_cloud_aerosol_per_sw_g_point) then
212      nspecsw = config%gas_optics_sw%spectral_def%ng
213    else
214      nspecsw = config%gas_optics_sw%spectral_def%nband
215    end if
216
217    if (config%do_cloud_aerosol_per_lw_g_point) then
218      nspeclw = config%gas_optics_lw%spectral_def%ng
219    else
220      nspeclw = config%gas_optics_lw%spectral_def%nband
221    end if
222
223    if (allocated(ao%wavelength_mono)) then
224      ! Monochromatic wavelengths also required
225      nmono = size(ao%wavelength_mono)
226    else
227      nmono = 0
228    end if
229
230    call ao%allocate(n_type_phobic, n_type_philic, nrh, nspeclw, nspecsw, nmono)
231
232    if (config%do_sw) then
233      call config%gas_optics_sw%spectral_def%calc_mapping(SolarReferenceTemperature, &
234           &  wavenumber, mapping, use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point))
235
236      ao%mass_ext_sw_phobic = matmul(mapping, mass_ext_phobic)
237      ao%ssa_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) &
238           &           / ao%mass_ext_sw_phobic
239      ao%g_sw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) &
240           &         / (ao%mass_ext_sw_phobic*ao%ssa_sw_phobic)
241
242      if (ao%use_hydrophilic) then
243        do jtype = 1,n_type_philic
244          ao%mass_ext_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype))
245          ao%ssa_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &
246               &                                        *ssa_philic(:,:,jtype)) &
247               &           / ao%mass_ext_sw_philic(:,:,jtype)
248          ao%g_sw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &
249               &                       *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) &
250               &         / (ao%mass_ext_sw_philic(:,:,jtype)*ao%ssa_sw_philic(:,:,jtype))
251        end do
252      end if
253    end if
254
255    if (config%do_lw) then
256      call config%gas_optics_lw%spectral_def%calc_mapping(TerrestrialReferenceTemperature, &
257           &  wavenumber, mapping, use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point))
258
259      ao%mass_ext_lw_phobic = matmul(mapping, mass_ext_phobic)
260      ao%ssa_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic) &
261           &           / ao%mass_ext_lw_phobic
262      ao%g_lw_phobic = matmul(mapping, mass_ext_phobic*ssa_phobic*g_phobic) &
263           &         / (ao%mass_ext_lw_phobic*ao%ssa_lw_phobic)
264
265      if (ao%use_hydrophilic) then
266        do jtype = 1,n_type_philic
267          ao%mass_ext_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype))
268          ao%ssa_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &
269               &                                        *ssa_philic(:,:,jtype)) &
270               &           / ao%mass_ext_lw_philic(:,:,jtype)
271          ao%g_lw_philic(:,:,jtype) = matmul(mapping, mass_ext_philic(:,:,jtype) &
272               &                       *ssa_philic(:,:,jtype)*g_philic(:,:,jtype)) &
273               &         / (ao%mass_ext_lw_philic(:,:,jtype)*ao%ssa_lw_philic(:,:,jtype))
274        end do
275      end if
276    end if
277
278    if (allocated(ao%wavelength_mono)) then
279      ! Monochromatic wavelengths also required
280      do jwl = 1,nmono
281        ! Wavelength (m) to wavenumber (cm-1)
282        wavenumber_target = 0.01_jprb / ao%wavelength_mono(jwl)
283        ! Find index to first interpolation point, and its weight
284        if (wavenumber_target <= wavenumber(1)) then
285          weight1 = 1.0_jprb
286          iwn = 1
287        else if (wavenumber_target >= wavenumber(nwn)) then
288          iwn = nwn-1
289          weight1 = 0.0_jprb
290        else
291          iwn = 1
292          do while (wavenumber(iwn+1) < wavenumber_target .and. iwn < nwn-1)
293            iwn = iwn + 1
294          end do
295          weight1 = (wavenumber(iwn+1)-wavenumber_target) &
296               &  / (wavenumber(iwn+1)-wavenumber(iwn))
297        end if
298        ! Linear interpolation
299        ao%mass_ext_mono_phobic(jwl,:) = weight1 * mass_ext_phobic(iwn,:) &
300             &             + (1.0_jprb - weight1)* mass_ext_phobic(iwn+1,:)
301        ao%ssa_mono_phobic(jwl,:)      = weight1 * ssa_phobic(iwn,:) &
302             &             + (1.0_jprb - weight1)* ssa_phobic(iwn+1,:)
303        ao%g_mono_phobic(jwl,:)        = weight1 * g_phobic(iwn,:) &
304             &             + (1.0_jprb - weight1)* g_phobic(iwn+1,:)
305        ao%lidar_ratio_mono_phobic(jwl,:) = weight1 * lidar_ratio_phobic(iwn,:) &
306             &                + (1.0_jprb - weight1)* lidar_ratio_phobic(iwn+1,:)
307        if (ao%use_hydrophilic) then
308          ao%mass_ext_mono_philic(jwl,:,:) = weight1 * mass_ext_philic(iwn,:,:) &
309               &               + (1.0_jprb - weight1)* mass_ext_philic(iwn+1,:,:)
310          ao%ssa_mono_philic(jwl,:,:)      = weight1 * ssa_philic(iwn,:,:) &
311               &               + (1.0_jprb - weight1)* ssa_philic(iwn+1,:,:)
312          ao%g_mono_philic(jwl,:,:)        = weight1 * g_philic(iwn,:,:) &
313               &               + (1.0_jprb - weight1)* g_philic(iwn+1,:,:)
314          ao%lidar_ratio_mono_philic(jwl,:,:) = weight1 * lidar_ratio_philic(iwn,:,:) &
315               &                  + (1.0_jprb - weight1)* lidar_ratio_philic(iwn+1,:,:)
316        end if
317      end do
318    end if
319
320    ! Deallocate memory local to this routine
321    deallocate(mass_ext_phobic)
322    deallocate(ssa_phobic)
323    deallocate(g_phobic)
324    deallocate(lidar_ratio_phobic)
325    if (ao%use_hydrophilic) then
326      deallocate(mass_ext_philic)
327      deallocate(ssa_philic)
328      deallocate(g_philic)
329      deallocate(lidar_ratio_philic)
330    end if
331
332    if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics',1,hook_handle)
333
334  end subroutine setup_general_aerosol_optics
335
336  !---------------------------------------------------------------------
337  ! Read file containing legacy-style band-wise aerosol optical
338  ! properties and average to the spectral intervals of the current
339  ! gas-optics scheme
340  subroutine setup_general_aerosol_optics_legacy(config, file_name)
341
342    use parkind1,                      only : jprb
343    use yomhook,                       only : lhook, dr_hook
344    use easy_netcdf,                   only : netcdf_file
345    use radiation_config,              only : config_type
346    use radiation_aerosol_optics_data, only : aerosol_optics_type
347    use radiation_spectral_definition, only : SolarReferenceTemperature, &
348         &                                    TerrestrialReferenceTemperature
349
350    type(config_type), intent(inout), target :: config
351
352    ! The NetCDF file containing the aerosol optics data
353    character(len=*), intent(in) :: file_name
354
355    ! Mapping matrix between optical properties at the wavenumbers in
356    ! the file, and spectral intervals used by the gas-optics scheme
357    real(jprb), allocatable :: mapping(:,:), mapping_transp(:,:)
358
359    ! Pointer to the aerosol optics coefficients for brevity of access
360    type(aerosol_optics_type), pointer :: ao
361
362    ! Local copy of aerosol optical properties in the spectral
363    ! intervals of the file, which is deallocated when it goes out of
364    ! scope
365    type(aerosol_optics_type) :: ao_legacy
366
367    integer :: jtype
368
369    real(jprb) :: hook_handle
370
371    if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics_legacy',0,hook_handle)
372    ao => config%aerosol_optics
373
374    ! Load file into a local structure
375    call ao_legacy%setup(file_name, iverbose=config%iverbosesetup)
376
377    ! Copy over scalars and coordinate variables
378    call ao%allocate(ao_legacy%n_type_phobic, ao_legacy%n_type_philic, ao_legacy%nrh, &
379         &           config%n_bands_lw, config%n_bands_sw, ao_legacy%n_mono_wl)
380    ao%description_phobic_str = ao_legacy%description_phobic_str
381    ao%description_philic_str = ao_legacy%description_philic_str
382    ao%rh_lower = ao_legacy%rh_lower
383
384    ! use_hydrophilic = ao_legacy%use_hydrophilic
385    ! ao%iclass = ao_legacy%iclass
386    ! ao%itype = ao_legacy%itype
387    ! ao%ntype = ao_legacy%ntype
388    ! ao%n_type_phobic = ao_legacy%n_type_phobic
389    ! ao%n_type_philic = ao_legacy%n_type_philic
390    ! ao%n_mono_wl = ao_legacy%n_mono_wl
391    ! ao%use_monochromatic = ao_legacy%use_monochromatic
392
393    if (config%do_sw) then
394      call config%gas_optics_sw%spectral_def%calc_mapping_from_wavenumber_bands(SolarReferenceTemperature, &
395           &  ao_legacy%wavenumber1_sw, ao_legacy%wavenumber2_sw, mapping_transp, &
396           &  use_bands=(.not. config%do_cloud_aerosol_per_sw_g_point))
397      if (allocated(mapping)) then
398        deallocate(mapping)
399      end if
400      allocate(mapping(config%n_bands_sw,ao_legacy%n_bands_sw))
401      mapping = transpose(mapping_transp)
402      ao%mass_ext_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic)
403      ao%ssa_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic*ao_legacy%ssa_sw_phobic) &
404           &           / ao%mass_ext_sw_phobic
405      ao%g_sw_phobic = matmul(mapping, ao_legacy%mass_ext_sw_phobic*ao_legacy%ssa_sw_phobic &
406           &                           *ao_legacy%g_sw_phobic) &
407           &         / (ao%mass_ext_sw_phobic*ao%ssa_sw_phobic)
408
409      if (ao%use_hydrophilic) then
410        do jtype = 1,ao%n_type_philic
411          ao%mass_ext_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype))
412          ao%ssa_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype) &
413               &                                        *ao_legacy%ssa_sw_philic(:,:,jtype)) &
414               &           / ao%mass_ext_sw_philic(:,:,jtype)
415          ao%g_sw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_sw_philic(:,:,jtype) &
416               &               *ao_legacy%ssa_sw_philic(:,:,jtype)*ao_legacy%g_sw_philic(:,:,jtype)) &
417               &         / (ao%mass_ext_sw_philic(:,:,jtype)*ao%ssa_sw_philic(:,:,jtype))
418        end do
419      end if
420    end if
421
422    if (config%do_lw) then
423      if (allocated(mapping_transp)) then
424        deallocate(mapping_transp)
425      end if
426      call config%gas_optics_lw%spectral_def%calc_mapping_from_wavenumber_bands(TerrestrialReferenceTemperature, &
427           &  ao_legacy%wavenumber1_lw, ao_legacy%wavenumber2_lw, mapping_transp, &
428           &  use_bands=(.not. config%do_cloud_aerosol_per_lw_g_point))
429      if (allocated(mapping)) then
430        deallocate(mapping)
431      end if
432      allocate(mapping(config%n_bands_lw,ao_legacy%n_bands_lw))
433      mapping = transpose(mapping_transp)
434      ao%mass_ext_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic)
435      ao%ssa_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic*ao_legacy%ssa_lw_phobic) &
436           &           / ao%mass_ext_lw_phobic
437      ao%g_lw_phobic = matmul(mapping, ao_legacy%mass_ext_lw_phobic*ao_legacy%ssa_lw_phobic &
438           &                           *ao_legacy%g_lw_phobic) &
439           &         / (ao%mass_ext_lw_phobic*ao%ssa_lw_phobic)
440
441      if (ao%use_hydrophilic) then
442        do jtype = 1,ao%n_type_philic
443          ao%mass_ext_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype))
444          ao%ssa_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype) &
445               &                                        *ao_legacy%ssa_lw_philic(:,:,jtype)) &
446               &           / ao%mass_ext_lw_philic(:,:,jtype)
447          ao%g_lw_philic(:,:,jtype) = matmul(mapping, ao_legacy%mass_ext_lw_philic(:,:,jtype) &
448               &               *ao_legacy%ssa_lw_philic(:,:,jtype)*ao_legacy%g_lw_philic(:,:,jtype)) &
449               &         / (ao%mass_ext_lw_philic(:,:,jtype)*ao%ssa_lw_philic(:,:,jtype))
450        end do
451      end if
452    end if
453
454    if (allocated(ao_legacy%wavelength_mono)) then
455      ao%wavelength_mono = ao_legacy%wavelength_mono
456      ao%mass_ext_mono_phobic = ao_legacy%mass_ext_mono_phobic
457      ao%ssa_mono_phobic = ao_legacy%ssa_mono_phobic
458      ao%g_mono_phobic = ao_legacy%g_mono_phobic
459      ao%lidar_ratio_mono_phobic = ao_legacy%lidar_ratio_mono_phobic
460      if (ao%use_hydrophilic) then
461        ao%mass_ext_mono_philic = ao_legacy%mass_ext_mono_philic
462        ao%ssa_mono_philic = ao_legacy%ssa_mono_philic
463        ao%g_mono_philic = ao_legacy%g_mono_philic
464        ao%lidar_ratio_mono_philic = ao_legacy%lidar_ratio_mono_philic
465      end if
466    end if
467
468    if (lhook) call dr_hook('radiation_aerosol_optics:setup_general_aerosol_optics_legacy',1,hook_handle)
469
470  end subroutine setup_general_aerosol_optics_legacy
471
472
473  !---------------------------------------------------------------------
474  ! Compute aerosol optical properties and add to existing gas optical
475  ! depth and scattering properties
476  subroutine add_aerosol_optics(nlev,istartcol,iendcol, &
477       &  config, thermodynamics, gas, aerosol, &
478       &  od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw)
479
480    use parkind1,                      only : jprb
481    use radiation_io,                  only : nulout, nulerr, radiation_abort
482    use yomhook,                       only : lhook, dr_hook
483    use radiation_config,              only : config_type
484    use radiation_thermodynamics,      only : thermodynamics_type
485    use radiation_gas,                 only : gas_type, IH2O, IMassMixingRatio
486    use radiation_aerosol,             only : aerosol_type
487    use radiation_constants,           only : AccelDueToGravity
488    use radiation_aerosol_optics_data, only : aerosol_optics_type, &
489         &  IAerosolClassUndefined,   IAerosolClassIgnored, &
490         &  IAerosolClassHydrophobic, IAerosolClassHydrophilic
491
492    integer, intent(in) :: nlev               ! number of model levels
493    integer, intent(in) :: istartcol, iendcol ! range of columns to process
494    type(config_type), intent(in), target :: config
495    type(thermodynamics_type),intent(in)  :: thermodynamics
496    type(gas_type),           intent(in)  :: gas
497    type(aerosol_type),       intent(in)  :: aerosol
498    ! Optical depth, single scattering albedo and asymmetry factor of
499    ! the atmosphere (gases on input, gases and aerosols on output)
500    ! for each g point. Note that longwave ssa and asymmetry and
501    ! shortwave asymmetry are all zero for gases, so are not yet
502    ! defined on input and are therefore intent(out).
503    real(jprb), dimension(config%n_g_lw,nlev,istartcol:iendcol), &
504         &   intent(inout) :: od_lw
505    real(jprb), dimension(config%n_g_lw_if_scattering,nlev,istartcol:iendcol), &
506         &   intent(out)   :: ssa_lw, g_lw
507    real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), &
508         &   intent(inout) :: od_sw, ssa_sw
509    real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), &
510         &   intent(out)   :: g_sw
511
512    ! Extinction optical depth, scattering optical depth and
513    ! asymmetry-times-scattering-optical-depth for all the aerosols at
514    ! a point in space for each spectral band of the shortwave and
515    ! longwave spectrum
516    real(jprb), dimension(config%n_bands_sw) &
517         & :: od_sw_aerosol, scat_sw_aerosol, scat_g_sw_aerosol, local_od_sw
518    real(jprb), dimension(config%n_bands_lw) :: od_lw_aerosol, local_od_lw
519    real(jprb), dimension(config%n_bands_lw_if_scattering) &
520         & :: scat_lw_aerosol, scat_g_lw_aerosol
521
522    real(jprb) :: h2o_mmr(istartcol:iendcol,nlev)
523
524    real(jprb) :: rh ! Relative humidity with respect to liquid water
525
526    ! Factor (kg m-2) to convert mixing ratio (kg kg-1) to mass in
527    ! path (kg m-2)
528    real(jprb) :: factor
529
530    ! Temporary extinction and scattering optical depths of aerosol
531    ! plus gas
532    real(jprb) :: local_od, local_scat
533
534    ! Loop indices for column, level, g point, band and aerosol type
535    integer :: jcol, jlev, jg, jtype, jband
536
537    ! Range of levels over which aerosols are present
538    integer :: istartlev, iendlev
539
540    ! Indices to spectral band and relative humidity look-up table
541    integer :: iband, irh
542
543    ! Short cut for ao%itype(jtype)
544    integer :: itype
545
546    ! Pointer to the aerosol optics coefficients for brevity of access
547    type(aerosol_optics_type), pointer :: ao
548
549    real(jprb) :: hook_handle
550
551    if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics',0,hook_handle)
552
553    if (aerosol%is_direct) then
554      ! Aerosol optical properties have been provided in each band
555      ! directly by the user
556      call add_aerosol_optics_direct(nlev,istartcol,iendcol, &
557           &  config, aerosol, &
558           &  od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw)
559    else
560      ! Aerosol mixing ratios have been provided
561
562      do jtype = 1,config%n_aerosol_types
563        if (config%aerosol_optics%iclass(jtype) == IAerosolClassUndefined) then
564          write(nulerr,'(a)') '*** Error: not all aerosol types are defined'
565          call radiation_abort()
566        end if
567      end do
568
569      if (config%iverbose >= 2) then
570        write(nulout,'(a)') 'Computing aerosol absorption/scattering properties'
571      end if
572
573      ao => config%aerosol_optics
574
575      istartlev = lbound(aerosol%mixing_ratio,2)
576      iendlev   = ubound(aerosol%mixing_ratio,2)
577
578      if (ubound(aerosol%mixing_ratio,3) /= config%n_aerosol_types) then
579        write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%mixing_ratio contains ', &
580             &  ubound(aerosol%mixing_ratio,3), ' aerosol types, expected ', &
581             &  config%n_aerosol_types
582        call radiation_abort()
583      end if
584
585      ! Set variables to zero that may not have been previously
586      g_sw(:,:,istartcol:iendcol) = 0.0_jprb
587      if (config%do_lw_aerosol_scattering) then
588        ssa_lw(:,:,istartcol:iendcol) = 0.0_jprb
589        g_lw(:,:,istartcol:iendcol)   = 0.0_jprb
590      end if
591
592      call gas%get(IH2O, IMassMixingRatio, h2o_mmr, istartcol=istartcol)
593
594      ! Loop over position
595      do jlev = istartlev,iendlev
596        do jcol = istartcol,iendcol
597          ! Compute relative humidity with respect to liquid
598          ! saturation and the index to the relative-humidity index of
599          ! hydrophilic-aerosol data
600          rh  = h2o_mmr(jcol,jlev) / thermodynamics%h2o_sat_liq(jcol,jlev)
601          irh = ao%calc_rh_index(rh)
602
603          factor = ( thermodynamics%pressure_hl(jcol,jlev+1) &
604               &    -thermodynamics%pressure_hl(jcol,jlev  )  ) &
605               &   / AccelDueToGravity 
606
607          ! Reset temporary arrays
608          od_sw_aerosol     = 0.0_jprb
609          scat_sw_aerosol   = 0.0_jprb
610          scat_g_sw_aerosol = 0.0_jprb
611          od_lw_aerosol     = 0.0_jprb
612          scat_lw_aerosol   = 0.0_jprb
613          scat_g_lw_aerosol = 0.0_jprb
614
615          do jtype = 1,config%n_aerosol_types
616
617            itype = ao%itype(jtype)
618
619            ! Add the optical depth, scattering optical depth and
620            ! scattering optical depth-weighted asymmetry factor for
621            ! this aerosol type to the total for all aerosols.  Note
622            ! that the following expressions are array-wise, the
623            ! dimension being spectral band.
624            if (ao%iclass(jtype) == IAerosolClassHydrophobic) then
625              do jband = 1,config%n_bands_sw
626                local_od_sw(jband) = factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
627                     &  * ao%mass_ext_sw_phobic(jband,itype)
628                od_sw_aerosol(jband) = od_sw_aerosol(jband) + local_od_sw(jband)
629                scat_sw_aerosol(jband) = scat_sw_aerosol(jband) &
630                     &  + local_od_sw(jband) * ao%ssa_sw_phobic(jband,itype)
631                scat_g_sw_aerosol(jband) = scat_g_sw_aerosol(jband) &
632                     &  + local_od_sw(jband) * ao%ssa_sw_phobic(jband,itype) &
633                     &  * ao%g_sw_phobic(jband,itype)
634              end do
635              if (config%do_lw_aerosol_scattering) then
636                local_od_lw = factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
637                     &  * ao%mass_ext_lw_phobic(:,itype)
638                od_lw_aerosol = od_lw_aerosol + local_od_lw
639                scat_lw_aerosol = scat_lw_aerosol &
640                     &  + local_od_lw * ao%ssa_lw_phobic(:,itype)
641                scat_g_lw_aerosol = scat_g_lw_aerosol &
642                     &  + local_od_lw * ao%ssa_lw_phobic(:,itype) &
643                     &  * ao%g_lw_phobic(:,itype)
644              else
645                ! If aerosol longwave scattering is not included then we
646                ! weight the optical depth by the single scattering
647                ! co-albedo
648                od_lw_aerosol = od_lw_aerosol &
649                     &  + factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
650                     &  * ao%mass_ext_lw_phobic(:,itype) &
651                     &  * (1.0_jprb - ao%ssa_lw_phobic(:,itype))
652              end if
653            else if (ao%iclass(jtype) == IAerosolClassHydrophilic) then
654              ! Hydrophilic aerosols require the look-up tables to
655              ! be indexed with irh
656              do jband = 1,config%n_bands_sw
657                local_od_sw(jband) = factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
658                     &  * ao%mass_ext_sw_philic(jband,irh,itype)
659                od_sw_aerosol(jband) = od_sw_aerosol(jband) + local_od_sw(jband)
660                scat_sw_aerosol(jband) = scat_sw_aerosol(jband) &
661                     &  + local_od_sw(jband) * ao%ssa_sw_philic(jband,irh,itype)
662                scat_g_sw_aerosol(jband) = scat_g_sw_aerosol(jband) &
663                     &  + local_od_sw(jband) * ao%ssa_sw_philic(jband,irh,itype) &
664                     &  * ao%g_sw_philic(jband,irh,itype)
665              end do
666              if (config%do_lw_aerosol_scattering) then
667                local_od_lw = factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
668                     &  * ao%mass_ext_lw_philic(:,irh,itype)
669                od_lw_aerosol = od_lw_aerosol + local_od_lw
670                scat_lw_aerosol = scat_lw_aerosol &
671                     &  + local_od_lw * ao%ssa_lw_philic(:,irh,itype)
672                scat_g_lw_aerosol = scat_g_lw_aerosol &
673                     &  + local_od_lw * ao%ssa_lw_philic(:,irh,itype) &
674                     &  * ao%g_lw_philic(:,irh,itype)
675              else
676                ! If aerosol longwave scattering is not included then we
677                ! weight the optical depth by the single scattering
678                ! co-albedo
679                od_lw_aerosol = od_lw_aerosol &
680                     &  + factor * aerosol%mixing_ratio(jcol,jlev,jtype) &
681                     &  * ao%mass_ext_lw_philic(:,irh,itype) &
682                     &  * (1.0_jprb - ao%ssa_lw_philic(:,irh,itype))
683              end if
684            end if
685            ! Implicitly, if ao%iclass(jtype) == IAerosolClassNone, then
686            ! no aerosol scattering properties are added
687
688          end do ! Loop over aerosol type
689
690          if (.not. config%do_sw_delta_scaling_with_gases) then
691            ! Delta-Eddington scaling on aerosol only.  Note that if
692            ! do_sw_delta_scaling_with_gases==.true. then the delta
693            ! scaling is done to the cloud-aerosol-gas mixture inside
694            ! the solver
695            call delta_eddington_extensive(od_sw_aerosol, scat_sw_aerosol, &
696                 &                         scat_g_sw_aerosol)
697          end if
698
699          ! Combine aerosol shortwave scattering properties with gas
700          ! properties (noting that any gas scattering will have an
701          ! asymmetry factor of zero)
702          do jg = 1,config%n_g_sw
703            iband = config%i_band_from_reordered_g_sw(jg)
704            local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband)
705            if (local_od > 0.0_jprb .and. od_sw_aerosol(iband) > 0.0_jprb) then
706              local_scat = ssa_sw(jg,jlev,jcol) * od_sw(jg,jlev,jcol) &
707                   &  + scat_sw_aerosol(iband)
708              ! Note that asymmetry_sw of gases is zero so the following
709              ! simply weights the aerosol asymmetry by the scattering
710              ! optical depth
711              if (local_scat > 0.0_jprb) then
712                g_sw(jg,jlev,jcol) = scat_g_sw_aerosol(iband) / local_scat
713              end if
714              ssa_sw(jg,jlev,jcol) = local_scat / local_od
715              od_sw (jg,jlev,jcol) = local_od
716            end if
717          end do
718
719          ! Combine aerosol longwave scattering properties with gas
720          ! properties, noting that in the longwave, gases do not
721          ! scatter at all
722          if (config%do_lw_aerosol_scattering) then
723
724            call delta_eddington_extensive(od_lw_aerosol, scat_lw_aerosol, &
725                 &                         scat_g_lw_aerosol) 
726
727            do jg = 1,config%n_g_lw
728              iband = config%i_band_from_reordered_g_lw(jg)
729              local_od = od_lw(jg,jlev,jcol) + od_lw_aerosol(iband)
730              if (local_od > 0.0_jprb .and. od_lw_aerosol(iband) > 0.0_jprb) then
731                ! All scattering is due to aerosols, therefore the
732                ! asymmetry factor is equal to the value for aerosols
733                if (scat_lw_aerosol(iband) > 0.0_jprb) then
734                  g_lw(jg,jlev,jcol) = scat_g_lw_aerosol(iband) &
735                       &  / scat_lw_aerosol(iband)
736                end if
737                ssa_lw(jg,jlev,jcol) = scat_lw_aerosol(iband) / local_od
738                od_lw (jg,jlev,jcol) = local_od
739              end if
740            end do
741          else
742            do jg = 1,config%n_g_lw
743              od_lw(jg,jlev,jcol) = od_lw(jg,jlev,jcol) &
744                   &  + od_lw_aerosol(config%i_band_from_reordered_g_lw(jg))
745            end do
746          end if
747
748        end do ! Loop over column
749      end do ! Loop over level
750
751    end if
752
753    if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics',1,hook_handle)
754
755  end subroutine add_aerosol_optics
756
757
758  !---------------------------------------------------------------------
759  ! Add precomputed optical properties to gas optical depth and
760  ! scattering properties
761  subroutine add_aerosol_optics_direct(nlev,istartcol,iendcol, &
762       &  config, aerosol, &
763       &  od_lw, ssa_lw, g_lw, od_sw, ssa_sw, g_sw)
764
765    use parkind1,                      only : jprb
766    use radiation_io,                  only : nulerr, radiation_abort
767    use yomhook,                       only : lhook, dr_hook
768    use radiation_config,              only : config_type
769    use radiation_aerosol,             only : aerosol_type
770
771    integer, intent(in) :: nlev               ! number of model levels
772    integer, intent(in) :: istartcol, iendcol ! range of columns to process
773    type(config_type), intent(in), target :: config
774    type(aerosol_type),       intent(in)  :: aerosol
775    ! Optical depth, single scattering albedo and asymmetry factor of
776    ! the atmosphere (gases on input, gases and aerosols on output)
777    ! for each g point. Note that longwave ssa and asymmetry and
778    ! shortwave asymmetry are all zero for gases, so are not yet
779    ! defined on input and are therefore intent(out).
780    real(jprb), dimension(config%n_g_lw,nlev,istartcol:iendcol), &
781         &   intent(inout) :: od_lw
782    real(jprb), dimension(config%n_g_lw_if_scattering,nlev,istartcol:iendcol), &
783         &   intent(out)   :: ssa_lw, g_lw
784    real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), &
785         &   intent(inout) :: od_sw, ssa_sw
786    real(jprb), dimension(config%n_g_sw,nlev,istartcol:iendcol), &
787         &   intent(out)   :: g_sw
788
789    ! Temporary extinction and scattering optical depths of aerosol
790    ! plus gas
791    real(jprb) :: local_od, local_scat
792
793    ! Extinction optical depth, scattering optical depth and
794    ! asymmetry-times-scattering-optical-depth for all the aerosols at
795    ! a point in space for each spectral band of the shortwave and
796    ! longwave spectrum
797    real(jprb), dimension(config%n_bands_sw,nlev) &
798         & :: od_sw_aerosol, scat_sw_aerosol, scat_g_sw_aerosol
799    real(jprb), dimension(config%n_bands_lw,nlev) :: od_lw_aerosol
800    real(jprb), dimension(config%n_bands_lw_if_scattering,nlev) &
801         & :: scat_lw_aerosol, scat_g_lw_aerosol
802
803    ! Loop indices for column, level, g point and band
804    integer :: jcol, jlev, jg, jb
805
806    ! Range of levels over which aerosols are present
807    integer :: istartlev, iendlev
808
809    ! Indices to spectral band
810    integer :: iband
811
812    real(jprb) :: hook_handle
813
814    if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics_direct',0,hook_handle)
815
816    if (config%do_sw) then
817      ! Check array dimensions
818      if (ubound(aerosol%od_sw,1) /= config%n_bands_sw) then
819        write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%od_sw contains ', &
820           &  ubound(aerosol%od_sw,1), ' band, expected ', &
821           &  config%n_bands_sw
822        call radiation_abort()
823      end if
824
825      istartlev = lbound(aerosol%od_sw,2)
826      iendlev   = ubound(aerosol%od_sw,2)
827
828      ! Set variables to zero that may not have been previously
829      g_sw(:,:,istartcol:iendcol) = 0.0_jprb
830
831      ! Loop over position
832      do jcol = istartcol,iendcol
833! Added for DWD (2020)
834!NEC$ forced_collapse
835        do jlev = istartlev,iendlev
836          do jb = 1,config%n_bands_sw
837            od_sw_aerosol(jb,jlev) = aerosol%od_sw(jb,jlev,jcol)
838            scat_sw_aerosol(jb,jlev) = aerosol%ssa_sw(jb,jlev,jcol) * od_sw_aerosol(jb,jlev)
839            scat_g_sw_aerosol(jb,jlev) = aerosol%g_sw(jb,jlev,jcol) * scat_sw_aerosol(jb,jlev)
840
841            if (.not. config%do_sw_delta_scaling_with_gases) then
842              ! Delta-Eddington scaling on aerosol only.  Note that if
843              ! do_sw_delta_scaling_with_gases==.true. then the delta
844              ! scaling is done to the cloud-aerosol-gas mixture
845              ! inside the solver
846              call delta_eddington_extensive(od_sw_aerosol(jb,jlev), scat_sw_aerosol(jb,jlev), &
847                   &                         scat_g_sw_aerosol(jb,jlev))
848            end if
849          end do
850        end do
851        ! Combine aerosol shortwave scattering properties with gas
852        ! properties (noting that any gas scattering will have an
853        ! asymmetry factor of zero)
854        do jlev = istartlev,iendlev
855          if (od_sw_aerosol(1,jlev) > 0.0_jprb) then
856            do jg = 1,config%n_g_sw
857              iband = config%i_band_from_reordered_g_sw(jg)
858              local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband,jlev)
859              local_scat = ssa_sw(jg,jlev,jcol) * od_sw(jg,jlev,jcol) &
860                   &  + scat_sw_aerosol(iband,jlev)
861              ! Note that asymmetry_sw of gases is zero so the following
862              ! simply weights the aerosol asymmetry by the scattering
863              ! optical depth
864              g_sw(jg,jlev,jcol) = scat_g_sw_aerosol(iband,jlev) / local_scat
865              local_od = od_sw(jg,jlev,jcol) + od_sw_aerosol(iband,jlev)
866              ssa_sw(jg,jlev,jcol) = local_scat / local_od
867              od_sw (jg,jlev,jcol) = local_od
868            end do
869          end if
870        end do
871      end do
872
873    end if
874
875    if (config%do_lw) then
876
877      if (ubound(aerosol%od_lw,1) /= config%n_bands_lw) then
878        write(nulerr,'(a,i0,a,i0)') '*** Error: aerosol%od_lw contains ', &
879           &  ubound(aerosol%od_lw,1), ' band, expected ', &
880           &  config%n_bands_lw
881        call radiation_abort()
882      end if
883
884      istartlev = lbound(aerosol%od_lw,2)
885      iendlev   = ubound(aerosol%od_lw,2)
886
887      if (config%do_lw_aerosol_scattering) then
888        ssa_lw(:,:,istartcol:iendcol) = 0.0_jprb
889        g_lw(:,:,istartcol:iendcol)   = 0.0_jprb
890 
891        ! Loop over position
892        do jcol = istartcol,iendcol
893! Added for DWD (2020)
894!NEC$ forced_collapse
895          do jlev = istartlev,iendlev
896            do jb = 1,config%n_bands_lw
897              od_lw_aerosol(jb,jlev) = aerosol%od_lw(jb,jlev,jcol)
898              scat_lw_aerosol(jb,jlev) = aerosol%ssa_lw(jb,jlev,jcol) * od_lw_aerosol(jb,jlev)
899              scat_g_lw_aerosol(jb,jlev) = aerosol%g_lw(jb,jlev,jcol) * scat_lw_aerosol(jb,jlev)
900           
901              call delta_eddington_extensive(od_lw_aerosol(jb,jlev), scat_lw_aerosol(jb,jlev), &
902                   &                         scat_g_lw_aerosol(jb,jlev))
903            end do
904          end do
905          do jlev = istartlev,iendlev
906            do jg = 1,config%n_g_lw
907              iband = config%i_band_from_reordered_g_lw(jg)
908              if (od_lw_aerosol(iband,jlev) > 0.0_jprb) then
909                ! All scattering is due to aerosols, therefore the
910                ! asymmetry factor is equal to the value for aerosols
911                if (scat_lw_aerosol(iband,jlev) > 0.0_jprb) then
912                  g_lw(jg,jlev,jcol) = scat_g_lw_aerosol(iband,jlev) &
913                       &  / scat_lw_aerosol(iband,jlev)
914                end if
915                local_od = od_lw(jg,jlev,jcol) + od_lw_aerosol(iband,jlev)
916                ssa_lw(jg,jlev,jcol) = scat_lw_aerosol(iband,jlev) / local_od
917                od_lw (jg,jlev,jcol) = local_od
918              end if
919            end do
920          end do
921        end do
922
923      else ! No longwave scattering
924
925        ! Loop over position
926        do jcol = istartcol,iendcol
927! Added for DWD (2020)
928!NEC$ forced_collapse
929          do jlev = istartlev,iendlev
930            ! If aerosol longwave scattering is not included then we
931            ! weight the optical depth by the single scattering
932            ! co-albedo
933            do jb = 1, config%n_bands_lw
934              od_lw_aerosol(jb,jlev) = aerosol%od_lw(jb,jlev,jcol) &
935                 &  * (1.0_jprb - aerosol%ssa_lw(jb,jlev,jcol))
936            end do
937          end do
938          do jlev = istartlev,iendlev
939            do jg = 1,config%n_g_lw
940              od_lw(jg,jlev,jcol) = od_lw(jg,jlev,jcol) &
941                   &  + od_lw_aerosol(config%i_band_from_reordered_g_lw(jg),jlev)
942            end do
943          end do
944        end do
945
946      end if
947    end if
948
949
950    if (lhook) call dr_hook('radiation_aerosol_optics:add_aerosol_optics_direct',1,hook_handle)
951
952  end subroutine add_aerosol_optics_direct
953 
954
955  !---------------------------------------------------------------------
956  ! Sometimes it is useful to specify aerosol in terms of its optical
957  ! depth at a particular wavelength.  This function returns the dry
958  ! mass-extinction coefficient, i.e. the extinction cross section per
959  ! unit mass, for aerosol of type "itype" at the specified wavelength
960  ! (m). For hydrophilic types, the value at the first relative
961  ! humidity bin is taken.
962  function dry_aerosol_mass_extinction(config, itype, wavelength)
963
964    use parkind1,                      only : jprb
965    use radiation_io,                  only : nulerr, radiation_abort
966    use radiation_config,              only : config_type
967    use radiation_aerosol_optics_data, only : aerosol_optics_type, &
968         &  IAerosolClassUndefined,   IAerosolClassIgnored, &
969         &  IAerosolClassHydrophobic, IAerosolClassHydrophilic
970
971    type(config_type), intent(in), target :: config
972
973    ! Aerosol type
974    integer, intent(in) :: itype
975
976    ! Wavelength (m)
977    real(jprb), intent(in) :: wavelength
978   
979    real(jprb) :: dry_aerosol_mass_extinction
980
981    ! Index to the monochromatic wavelength requested
982    integer :: imono
983
984    ! Pointer to the aerosol optics coefficients for brevity of access
985    type(aerosol_optics_type), pointer :: ao
986
987    ao => config%aerosol_optics
988
989    imono = minloc(abs(wavelength - ao%wavelength_mono), 1)
990
991    if (abs(wavelength - ao%wavelength_mono(imono))/wavelength > 0.01_jprb) then
992      write(nulerr,'(a,e8.4,a)') '*** Error: requested wavelength ', &
993           &  wavelength, ' not within 1% of stored wavelengths'
994      call radiation_abort()
995     end if
996
997    if (ao%iclass(itype) == IAerosolClassHydrophobic) then
998      dry_aerosol_mass_extinction = ao%mass_ext_mono_phobic(imono,ao%itype(itype))
999    else if (ao%iclass(itype) == IAerosolClassHydrophilic) then
1000      ! Take the value at the first relative-humidity bin for the
1001      ! "dry" aerosol value
1002      dry_aerosol_mass_extinction = ao%mass_ext_mono_philic(imono,1,ao%itype(itype))
1003    else
1004      dry_aerosol_mass_extinction = 0.0_jprb
1005    end if
1006
1007  end function dry_aerosol_mass_extinction
1008
1009
1010  !---------------------------------------------------------------------
1011  ! Compute aerosol extinction coefficient at a particular wavelength
1012  ! and a single height - this is useful for visibility diagnostics
1013  subroutine aerosol_extinction(ncol,istartcol,iendcol, &
1014       &  config, wavelength, mixing_ratio, relative_humidity, extinction)
1015
1016    use parkind1,                      only : jprb
1017    use yomhook,                       only : lhook, dr_hook
1018    use radiation_io,                  only : nulerr, radiation_abort
1019    use radiation_config,              only : config_type
1020    use radiation_aerosol_optics_data, only : aerosol_optics_type, &
1021         &  IAerosolClassUndefined,   IAerosolClassIgnored, &
1022         &  IAerosolClassHydrophobic, IAerosolClassHydrophilic
1023
1024    integer, intent(in) :: ncol               ! number of columns
1025    integer, intent(in) :: istartcol, iendcol ! range of columns to process
1026    type(config_type), intent(in), target :: config
1027    real(jprb), intent(in)  :: wavelength ! Requested wavelength (m)
1028    real(jprb), intent(in)  :: mixing_ratio(ncol,config%n_aerosol_types)
1029    real(jprb), intent(in)  :: relative_humidity(ncol)
1030    real(jprb), intent(out) :: extinction(ncol)
1031
1032    ! Local aerosol extinction
1033    real(jprb) :: ext
1034
1035    ! Index to the monochromatic wavelength requested
1036    integer :: imono
1037
1038    ! Pointer to the aerosol optics coefficients for brevity of access
1039    type(aerosol_optics_type), pointer :: ao
1040   
1041    ! Loop indices for column and aerosol type
1042    integer :: jcol, jtype
1043
1044    ! Relative humidity index
1045    integer :: irh
1046
1047    real(jprb) :: hook_handle
1048
1049    if (lhook) call dr_hook('radiation_aerosol_optics:aerosol_extinction',0,hook_handle)
1050
1051    do jtype = 1,config%n_aerosol_types
1052      if (config%aerosol_optics%iclass(jtype) == IAerosolClassUndefined) then
1053        write(nulerr,'(a)') '*** Error: not all aerosol types are defined'
1054        call radiation_abort()
1055      end if
1056    end do
1057
1058    ao => config%aerosol_optics
1059
1060    imono = minloc(abs(wavelength - ao%wavelength_mono), 1)
1061
1062    if (abs(wavelength - ao%wavelength_mono(imono))/wavelength > 0.01_jprb) then
1063      write(nulerr,'(a,e8.4,a)') '*** Error: requested wavelength ', &
1064           &  wavelength, ' not within 1% of stored wavelengths'
1065      call radiation_abort()
1066     end if
1067
1068    ! Loop over position
1069    do jcol = istartcol,iendcol
1070      ext = 0.0_jprb
1071      ! Get relative-humidity index
1072      irh = ao%calc_rh_index(relative_humidity(jcol))
1073      ! Add extinction coefficients from each aerosol type
1074      do jtype = 1,config%n_aerosol_types
1075        if (ao%iclass(jtype) == IAerosolClassHydrophobic) then
1076          ext = ext + mixing_ratio(jcol,jtype) &
1077               &    * ao%mass_ext_mono_phobic(imono,ao%itype(jtype))
1078        else if (ao%iclass(jtype) == IAerosolClassHydrophilic) then
1079          ext = ext + mixing_ratio(jcol,jtype) &
1080               &    * ao%mass_ext_mono_philic(imono,irh,ao%itype(jtype))
1081        end if
1082      end do
1083
1084      extinction(jcol) = ext
1085    end do
1086
1087    if (lhook) call dr_hook('radiation_aerosol_optics:aerosol_extinction',1,hook_handle)
1088
1089  end subroutine aerosol_extinction
1090
1091end module radiation_aerosol_optics
Note: See TracBrowser for help on using the repository browser.