mam  v1.0
A Modal Aerosol Model
optics_util.F90
Go to the documentation of this file.
1! Copyright (C) 2021 National Center for Atmospheric Research
2! SPDX-License-Identifier: Apache-2.0
3!
4!> \file
5!> The mam_optics_constants module
6
7!> Constants used to calculate MAM optical properties
9
10 use ai_optics, only : optics_t
11 use ai_optics_absorption_optical_depth, &
12 only : optics_absorption_optical_depth_t
13 use ai_optics_asymmetry_factor, only : optics_asymmetry_factor_t
14 use ai_optics_extinction_optical_depth, &
15 only : optics_extinction_optical_depth_t
16 use ai_optics_forward_scattered_fraction, &
17 only : optics_forward_scattered_fraction_t
18 use ai_optics_single_scatter_albedo, &
19 only : optics_single_scatter_albedo_t
20 use musica_constants, only : musica_dk
21
22 implicit none
23 private
24
26
27contains
28
29!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
30
31 !> Returns an optics_t object for a given property
32 function create_optics( property, native_shortwave_grid, &
33 native_longwave_grid, output_grid, interpolation_strategy ) &
34 result( new_optics )
35
36 use ai_wavelength_grid, only : wavelength_grid_t, kwavenumber, &
37 kcentimeter
38 use musica_assert, only : die_msg
39 use musica_interpolator, only : interpolation_strategy_i
40 use musica_property, only : property_t
41 use musica_string, only : string_t
42
43 class(optics_t), pointer :: new_optics
44 class(property_t), intent(in) :: property
45 type(wavelength_grid_t), intent(in) :: native_shortwave_grid
46 type(wavelength_grid_t), intent(in) :: native_longwave_grid
47 type(wavelength_grid_t), intent(in) :: output_grid
48 procedure(interpolation_strategy_i), optional :: interpolation_strategy
49
50 type(wavelength_grid_t) :: native_grid
51 type(string_t) :: property_name
52
53 property_name = property%name( )
54 if( property_name .eq. "layer extinction optical depth" ) then
55 new_optics => optics_extinction_optical_depth_t( native_shortwave_grid, &
56 output_grid, &
57 interpolation_strategy = interpolation_strategy )
58 else if( property_name .eq. "layer single-scatter albedo" ) then
59 new_optics => optics_single_scatter_albedo_t( native_shortwave_grid, &
60 output_grid, &
61 interpolation_strategy = interpolation_strategy )
62 else if( property_name .eq. "asymmetry factor" ) then
63 new_optics => optics_asymmetry_factor_t( native_shortwave_grid, &
64 output_grid, &
65 interpolation_strategy = interpolation_strategy )
66 else if( property_name .eq. "forward scattered fraction" ) then
67 new_optics => optics_forward_scattered_fraction_t( &
68 native_shortwave_grid, &
69 output_grid, &
70 interpolation_strategy = interpolation_strategy )
71 else if( property_name .eq. "layer absorption optical depth" ) then
72 new_optics => optics_absorption_optical_depth_t( native_longwave_grid, &
73 output_grid, &
74 interpolation_strategy = interpolation_strategy )
75 else
76 call die_msg( 769442313, "Unsupported optical property: '"// &
77 property_name//"'" )
78 end if
79
80 end function create_optics
81
82!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
83
84 !> Adds extinction optical depths to optical property values
86 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
87 optics )
88
89 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
90 real(kind=musica_dk), intent(in) :: single_scatter_albedo(:)
91 real(kind=musica_dk), intent(in) :: asymmetry_factor(:)
92 class(optics_extinction_optical_depth_t), intent(inout) :: optics
93
94 call optics%add_values( extinction_optical_depth )
95
97
98!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
99
100 !> Adds single-scatter albedo to optical property values
102 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
103 optics )
104
105 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
106 real(kind=musica_dk), intent(in) :: single_scatter_albedo(:)
107 real(kind=musica_dk), intent(in) :: asymmetry_factor(:)
108 class(optics_single_scatter_albedo_t), intent(inout) :: optics
109
110 !> \todo the parameter calculated from the lookup-table is called
111 !! single-scatter albedo, as is the value returned to radiation
112 !! after multiplying by extinction. Which is correct?
113 call optics%add_values( extinction_optical_depth(:) &
114 * single_scatter_albedo(:) )
115
117
118!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
119
120 !> Adds asymmetry factor to optical property values
122 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
123 optics )
124
125 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
126 real(kind=musica_dk), intent(in) :: single_scatter_albedo(:)
127 real(kind=musica_dk), intent(in) :: asymmetry_factor(:)
128 class(optics_asymmetry_factor_t), intent(inout) :: optics
129
130 !> \todo the parameter calculated from the lookup-table is called
131 !! asymmetry factor, as is the value returned to radiation
132 !! after multiplying by extinction and single-scatter albedo.
133 !! Which is correct?
134 call optics%add_values( extinction_optical_depth(:) &
135 * single_scatter_albedo(:) &
136 * asymmetry_factor(:) )
137
138 end subroutine add_shortwave_asymmetry_factor
139
140!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
141
142 !> Adds forward scattered fraction to optical property values
144 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
145 optics )
146
147 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
148 real(kind=musica_dk), intent(in) :: single_scatter_albedo(:)
149 real(kind=musica_dk), intent(in) :: asymmetry_factor(:)
150 class(optics_forward_scattered_fraction_t), intent(inout) :: optics
151
152 call optics%add_values( extinction_optical_depth(:) &
153 * single_scatter_albedo(:) &
154 * asymmetry_factor(:) &
155 * asymmetry_factor(:) )
156
158
159!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
160
161 !> Adds to shortwave property values
163 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
164 optics )
165
166 use musica_assert, only : die_msg
167
168 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
169 real(kind=musica_dk), intent(in) :: single_scatter_albedo(:)
170 real(kind=musica_dk), intent(in) :: asymmetry_factor(:)
171 class(optics_t), intent(inout) :: optics
172
173 select type( optics )
174 class is( optics_extinction_optical_depth_t )
176 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
177 optics )
178 class is( optics_single_scatter_albedo_t )
180 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
181 optics )
182 class is( optics_asymmetry_factor_t )
184 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
185 optics )
186 class is( optics_forward_scattered_fraction_t )
188 extinction_optical_depth, single_scatter_albedo, asymmetry_factor, &
189 optics )
190 class default
191 call die_msg( 628273876, "Unsupported MAM shortwave optical property" )
192 end select
193
194 end subroutine add_shortwave_property
195
196!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
197
198 !> Adds absorption optical depth to optical property values
200 extinction_optical_depth, optics )
201
202 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
203 class(optics_absorption_optical_depth_t), intent(inout) :: optics
204
205 call optics%add_values( extinction_optical_depth )
206
208
209!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
210
211 !> Adds to longwave property values
212 subroutine add_longwave_property( extinction_optical_depth, optics )
213
214 use musica_assert, only : die_msg
215
216 real(kind=musica_dk), intent(in) :: extinction_optical_depth(:)
217 class(optics_t), intent(inout) :: optics
218
219 select type( optics )
220 class is( optics_absorption_optical_depth_t )
221 call add_longwave_absorption_optical_depth( extinction_optical_depth, &
222 optics )
223 class default
224 call die_msg(405116022, "Unsupported MAM longwave optical property" )
225 end select
226
227 end subroutine add_longwave_property
228
229!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
230
231end module mam_optics_util
Constants used to calculate MAM optical properties.
Definition: optics_util.F90:8
subroutine add_shortwave_asymmetry_factor(extinction_optical_depth, single_scatter_albedo, asymmetry_factor, optics)
Adds asymmetry factor to optical property values.
subroutine add_shortwave_extinction_optical_depth(extinction_optical_depth, single_scatter_albedo, asymmetry_factor, optics)
Adds extinction optical depths to optical property values.
Definition: optics_util.F90:88
subroutine add_longwave_absorption_optical_depth(extinction_optical_depth, optics)
Adds absorption optical depth to optical property values.
subroutine add_shortwave_forward_scattered_fraction(extinction_optical_depth, single_scatter_albedo, asymmetry_factor, optics)
Adds forward scattered fraction to optical property values.
subroutine add_shortwave_single_scatter_albedo(extinction_optical_depth, single_scatter_albedo, asymmetry_factor, optics)
Adds single-scatter albedo to optical property values.
class(optics_t) function, pointer, public create_optics(property, native_shortwave_grid, native_longwave_grid, output_grid, interpolation_strategy)
Returns an optics_t object for a given property.
Definition: optics_util.F90:35
subroutine, public add_shortwave_property(extinction_optical_depth, single_scatter_albedo, asymmetry_factor, optics)
Adds to shortwave property values.
subroutine, public add_longwave_property(extinction_optical_depth, optics)
Adds to longwave property values.
real(kind=musica_dk), dimension(:,:), intent(out), optional asymmetry_factor