mam  v1.0
A Modal Aerosol Model
core.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_core module
6
7!> The core_t type and related functions
8module mam_core
9
10 use ai_aerosol, only : aerosol_t
11 use ai_aerosol_state, only : aerosol_state_t
12 use ai_wavelength_grid, only : wavelength_grid_t
13 use mam_mode, only : mode_t, mode_state_t
14
15 implicit none
16 private
17
18 public :: core_t, state_t
19
20 !> The Modal Aerosol Model core
21 type, extends(aerosol_t) :: core_t
22 private
23 type(wavelength_grid_t) :: shortwave_grid_
24 type(wavelength_grid_t) :: longwave_grid_
25 type(mode_t), allocatable :: modes_(:)
26 contains
27 procedure :: new_state
28 procedure :: new_optics
29 procedure, private :: shortwave_optics_scalar
30 procedure, private :: shortwave_optics_array
31 procedure, private :: longwave_optics_scalar
32 procedure, private :: longwave_optics_array
33 procedure :: print_state
34 end type core_t
35
36 interface core_t
37 module procedure :: constructor
38 end interface
39
40 !> Modal aerosol state
41 type, extends(aerosol_state_t) :: state_t
42 private
43 type(mode_state_t), allocatable :: mode_states_(:)
44 contains
45 procedure :: raw_size
46 procedure :: load_state
47 procedure :: dump_state
48 procedure :: randomize
49 end type state_t
50
51contains
52
53!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
54
55 !> Constructor of the MAM core
56 function constructor( config ) result( new_core )
57
58 use mam_species, only : species_t
59 use musica_assert, only : assert_msg
60 use musica_config, only : config_t
61 use musica_iterator, only : iterator_t
62 use musica_string, only : string_t
63
64 type(core_t), pointer :: new_core
65 class(config_t), intent(inout) :: config
66
67 character(len=*), parameter :: my_name = "MAM core_t constructor"
68 type(config_t) :: modes, mode
69 type(string_t) :: file_name
70 class(iterator_t), pointer :: iter
71 integer :: i_mode
72 type(wavelength_grid_t) :: mode_grid
73
74 allocate( new_core )
75 call config%get( "modes", modes, my_name )
76 allocate( new_core%modes_( modes%number_of_children( ) ) )
77 iter => modes%get_iterator( )
78 i_mode = 1
79 do while( iter%next( ) )
80 call modes%get( iter, mode, my_name )
81 new_core%modes_( i_mode ) = mode_t( mode )
82 if( i_mode .eq. 1 ) then
83 new_core%shortwave_grid_ = new_core%modes_( i_mode )%shortwave_grid( )
84 new_core%longwave_grid_ = new_core%modes_( i_mode )%longwave_grid( )
85 else
86 call assert_msg( 253107950, new_core%shortwave_grid_ .eq. &
87 new_core%modes_( i_mode )%shortwave_grid( ), &
88 "MAM is not currently set up to use different "// &
89 "shortwave wavelength grids for individual modes." )
90 call assert_msg( 525506880, new_core%longwave_grid_ .eq. &
91 new_core%modes_( i_mode )%longwave_grid( ), &
92 "MAM is not currently set up to use different "// &
93 "longwave wavelength grids for individual modes." )
94 end if
95 i_mode = i_mode + 1
96 end do
97 deallocate( iter )
98
99 end function constructor
100
101!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
102
103 !> Creates a new state object for the modal aerosol
104 function new_state( this )
105
106 class(aerosol_state_t), pointer :: new_state
107 class(core_t), intent(in) :: this
108
109 class(aerosol_state_t), pointer :: mode_state
110 integer :: i_mode
111
112 allocate( state_t :: new_state )
113 select type( new_state )
114 type is( state_t )
115 allocate( new_state%mode_states_( size( this%modes_ ) ) )
116 do i_mode = 1, size( this%modes_ )
117 mode_state => this%modes_( i_mode )%new_state( )
118 select type( mode_state )
119 class is( mode_state_t )
120 new_state%mode_states_( i_mode ) = mode_state
121 end select
122 deallocate( mode_state )
123 end do
124 end select
125
126 end function new_state
127
128!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
129
130 !> Creates an optics_t object for a given property
131 function new_optics( this, property, output_grid, interpolation_strategy )
132
133 use ai_optics, only : optics_t
134 use ai_wavelength_grid, only : wavelength_grid_t
136 use musica_interpolator, only : interpolation_strategy_i
137 use musica_property, only : property_t
138
139 class(optics_t), pointer :: new_optics
140 class(core_t), intent(in) :: this
141 class(property_t), intent(in) :: property
142 type(wavelength_grid_t), intent(in) :: output_grid
143 procedure(interpolation_strategy_i), optional :: interpolation_strategy
144
145 new_optics => &
146 create_optics( property, this%shortwave_grid_, this%longwave_grid_, &
147 output_grid, interpolation_strategy )
148
149 end function new_optics
150
151!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
152
153 !> Returns shortwave optical properties
154 subroutine shortwave_optics_scalar( this, environmental_state, &
155 aerosol_state, optics )
156
157 use ai_environmental_state, only : environmental_state_t
158 use ai_optics, only : optics_t, optics_ptr
159
160 class(core_t), intent(in) :: this
161 class(environmental_state_t), intent(in) :: environmental_state
162 class(aerosol_state_t), intent(in) :: aerosol_state
163 class(optics_t), target, intent(inout) :: optics
164
165 type(optics_ptr) :: optics_set(1)
166
167 optics_set(1)%ptr_ => optics
168 call this%shortwave_optics( environmental_state, aerosol_state, &
169 optics_set )
170 nullify( optics_set(1)%ptr_ )
171
172 end subroutine shortwave_optics_scalar
173
174!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
175
176 !> Returns shortwave optical properties
177 subroutine shortwave_optics_array( this, environmental_state, &
178 aerosol_state, optics )
179
180 use ai_environmental_state, only : environmental_state_t
181 use ai_optics, only : optics_ptr
182 use musica_assert, only : die_msg
183
184 class(core_t), intent(in) :: this
185 class(environmental_state_t), intent(in) :: environmental_state
186 class(aerosol_state_t), intent(in) :: aerosol_state
187 type(optics_ptr), intent(inout) :: optics(:)
188
189 integer :: i_mode, i_prop
190
191 select type( aerosol_state )
192 class is( state_t )
193 do i_prop = 1, size( optics )
194 call optics( i_prop )%ptr_%reset_values( )
195 end do
196 do i_mode = 1, size( this%modes_ )
197 associate( mode_state => aerosol_state%mode_states_( i_mode ), &
198 mode => this%modes_( i_mode ) )
199 call mode%add_shortwave_optics( environmental_state, &
200 mode_state, &
201 optics )
202 end associate
203 end do
204 class default
205 call die_msg( 822004975, "Invalid state passed to MAM shortwave "// &
206 "optics calculator" )
207 end select
208
209 end subroutine shortwave_optics_array
210
211!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
212
213 !> Returns longwave optical properties
214 subroutine longwave_optics_scalar( this, environmental_state, &
215 aerosol_state, optics )
216
217 use ai_environmental_state, only : environmental_state_t
218 use ai_optics, only : optics_t, optics_ptr
219
220 class(core_t), intent(in) :: this
221 class(environmental_state_t), intent(in) :: environmental_state
222 class(aerosol_state_t), intent(in) :: aerosol_state
223 class(optics_t), target, intent(inout) :: optics
224
225 type(optics_ptr) :: optics_set(1)
226
227 optics_set(1)%ptr_ => optics
228 call this%longwave_optics( environmental_state, aerosol_state, &
229 optics_set )
230 nullify( optics_set(1)%ptr_ )
231
232 end subroutine longwave_optics_scalar
233
234!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
235
236 !> Returns longwave optical properties
237 subroutine longwave_optics_array( this, environmental_state, &
238 aerosol_state, optics )
239
240 use ai_environmental_state, only : environmental_state_t
241 use ai_optics, only : optics_ptr
242 use musica_assert, only : die_msg
243
244 class(core_t), intent(in) :: this
245 class(environmental_state_t), intent(in) :: environmental_state
246 class(aerosol_state_t), intent(in) :: aerosol_state
247 type(optics_ptr), intent(inout) :: optics(:)
248
249 integer :: i_mode, i_prop
250
251 select type( aerosol_state )
252 class is( state_t )
253 do i_prop = 1, size( optics )
254 call optics( i_prop )%ptr_%reset_values( )
255 end do
256 do i_mode = 1, size( this%modes_ )
257 associate( mode_state => aerosol_state%mode_states_( i_mode ), &
258 mode => this%modes_( i_mode ) )
259 call mode%add_longwave_optics( environmental_state, &
260 mode_state, &
261 optics )
262 end associate
263 end do
264 class default
265 call die_msg( 466350541, "Invalid state passed to MAM longwave "// &
266 "optics calculator" )
267 end select
268
269 end subroutine longwave_optics_array
270
271!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
272
273 !> Ouptuts the current aerosol state
274 subroutine print_state( this, aerosol_state, io_unit )
275
276 use musica_assert, only : die_msg
277 use musica_string, only : to_char
278
279 !> MAM core
280 class(core_t), intent(in) :: this
281 !> MAM aerosol state
282 class(aerosol_state_t), intent(in) :: aerosol_state
283 !> Optional output unit (defaults to 6)
284 integer, optional, intent(in) :: io_unit
285
286 integer :: lunit, i_mode
287
288 lunit = 6
289 if( present( io_unit ) ) lunit = io_unit
290 select type( aerosol_state )
291 class is( state_t )
292 write(lunit,*) "**** MAM Aerosol State ****"
293 if( .not. allocated( this%modes_ ) ) then
294 write(lunit,*) "--- Uninitialized MAM core ---"
295 write(lunit,*) "**** End MAM Aerosol State ****"
296 return
297 end if
298 if( .not. allocated( aerosol_state%mode_states_ ) ) then
299 write(lunit,*) "--- Uninitialized MAM state ---"
300 write(lunit,*) "**** End MAM Aerosol State ****"
301 return
302 end if
303 do i_mode = 1, size( this%modes_ )
304 write(lunit,*) "*** Mode "//trim( to_char( i_mode ) )//" ***"
305 call this%modes_( i_mode )%print_state( &
306 aerosol_state%mode_states_( i_mode ), io_unit )
307 end do
308 write(lunit,*) "**** End MAM Aerosol State ****"
309 class default
310 call die_msg( 970908819, "Invalid state passed to MAM aerosol" )
311 end select
312
313 end subroutine print_state
314
315!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
316!!
317!! state_t functions
318!!
319!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
320
321 !> Returns the number of doubles needed to hold the raw MAM aerosol state
322 integer function raw_size( this )
323
324 class(state_t), intent(in) :: this
325
326 integer :: i_mode
327
328 raw_size = 0
329 do i_mode = 1, size( this%mode_states_ )
330 raw_size = raw_size + this%mode_states_( i_mode )%raw_size( )
331 end do
332
333 end function raw_size
334
335!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
336
337 !> Loads raw MAM state data to the state_t object
338 subroutine load_state( this, raw_state, index )
339
340 use musica_constants, only : musica_dk
341
342 class(state_t), intent(inout) :: this
343 real(kind=musica_dk), intent(in) :: raw_state(:)
344 integer, optional, intent(inout) :: index
345
346 integer :: i_mode, lindex
347
348 lindex = 1
349 if( present( index ) ) lindex = index
350 do i_mode = 1, size( this%mode_states_ )
351 call this%mode_states_( i_mode )%load_state( raw_state, lindex )
352 end do
353 if( present( index ) ) index = lindex
354
355 end subroutine load_state
356
357!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
358
359 !> Dumps the raw MAM state data to an array
360 subroutine dump_state( this, raw_state, index )
361
362 use musica_constants, only : musica_dk
363
364 class(state_t), intent(in) :: this
365 real(kind=musica_dk), intent(inout) :: raw_state(:)
366 integer, optional, intent(inout) :: index
367
368 integer :: i_mode, lindex
369
370 lindex = 1
371 if( present( index ) ) lindex = index
372 do i_mode = 1, size( this%mode_states_ )
373 call this%mode_states_( i_mode )%dump_state( raw_state, lindex )
374 end do
375 if( present( index ) ) index = lindex
376
377 end subroutine dump_state
378
379!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
380
381 !> Sets the MAM state to a random, but reasonable, state. For testing only.
382 subroutine randomize( this )
383
384 use musica_assert, only : assert_msg
385
386 class(state_t), intent(inout) :: this
387
388 integer :: i_mode
389
390 call assert_msg( 652075662, allocated( this%mode_states_ ), &
391 "Trying to randomize uninitialized MAM state" )
392 do i_mode = 1, size( this%mode_states_ )
393 call this%mode_states_( i_mode )%randomize( )
394 end do
395
396 end subroutine randomize
397
398!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
399
400end module mam_core
The core_t type and related functions.
Definition: core.F90:8
subroutine print_state(this, aerosol_state, io_unit)
Ouptuts the current aerosol state.
Definition: core.F90:275
subroutine randomize(this)
Sets the MAM state to a random, but reasonable, state. For testing only.
Definition: core.F90:383
subroutine load_state(this, raw_state, index)
Loads raw MAM state data to the state_t object.
Definition: core.F90:339
subroutine shortwave_optics_scalar(this, environmental_state, aerosol_state, optics)
Returns shortwave optical properties.
Definition: core.F90:156
type(core_t) function, pointer constructor(config)
Constructor of the MAM core.
Definition: core.F90:57
class(optics_t) function, pointer new_optics(this, property, output_grid, interpolation_strategy)
Creates an optics_t object for a given property.
Definition: core.F90:132
subroutine longwave_optics_scalar(this, environmental_state, aerosol_state, optics)
Returns longwave optical properties.
Definition: core.F90:216
integer function raw_size(this)
Returns the number of doubles needed to hold the raw MAM aerosol state.
Definition: core.F90:323
subroutine longwave_optics_array(this, environmental_state, aerosol_state, optics)
Returns longwave optical properties.
Definition: core.F90:239
class(aerosol_state_t) function, pointer new_state(this)
Creates a new state object for the modal aerosol.
Definition: core.F90:105
subroutine shortwave_optics_array(this, environmental_state, aerosol_state, optics)
Returns shortwave optical properties.
Definition: core.F90:179
subroutine dump_state(this, raw_state, index)
Dumps the raw MAM state data to an array.
Definition: core.F90:361
The mode_t type and related functions.
Definition: mode.F90:8
Constants used to calculate MAM optical properties.
Definition: optics_util.F90:8
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
The species_t type and related functions.
Definition: species.F90:8
The Modal Aerosol Model core.
Definition: core.F90:21
Modal aerosol state.
Definition: core.F90:41
Aerosol mode state.
Definition: mode.F90:66
An aerosol mode.
Definition: mode.F90:22