mam  v1.0
A Modal Aerosol Model
Data Types | Functions/Subroutines
mam_mode Module Reference

The mode_t type and related functions. More...

Data Types

type  mode_state_t
 Aerosol mode state. More...
 
interface  mode_t
 An aerosol mode. More...
 

Functions/Subroutines

subroutine add_longwave_optics (this, environmental_state, mode_state, optics)
 Calculates longwave optical properties for a given state and adds them to the set of optical properties. More...
 
subroutine add_shortwave_optics (this, environmental_state, mode_state, optics)
 Calculates shortwave optical properties for a given state and adds them to the set of optical properties. More...
 
pure real(kind=musica_dk) function, dimension(number_of_bands) asymmetry_factor (this, mode_state, number_of_bands, number_of_coefficients, coefficients, size_function)
 Calculates the asymmetry parameter [unitless ratio] for a given Chebyshev function. More...
 
type(mode_t) function constructor (config)
 Constructor of mode_t objects. More...
 
subroutine dump_state (this, raw_state, index)
 Dumps the aerosol state into the raw state array. More...
 
elemental real(kind=musica_dk) function geometric_mean_diameter_of_number_distribution__m (this, mode_state)
 Returns the geometric mean diameter [m] of the number distribution of the mode. More...
 
elemental real(kind=musica_dk) function geometric_standard_deviation_of_number_distribution (this, mode_state)
 Returns the geometric standard deviation of the number distribution of the mode. More...
 
subroutine load_state (this, raw_state, index)
 Loads raw mode state data to the mode_state_t object. More...
 
type(wavelength_grid_t) function longwave_grid (this)
 Returns the native longwave optics wavelength grid. More...
 
subroutine longwave_optics_array (this, environmental_state, aerosol_state, optics)
 Returns longwave optical properties. More...
 
subroutine longwave_optics_scalar (this, environmental_state, aerosol_state, optics)
 Returns longwave optical properties. More...
 
pure complex(kind=musica_dk) function, dimension(number_of_bands) net_longwave_refractive_index (this, mode_state, number_of_bands)
 Calculates the net refractive index for the mode based on weighted contributions from constituent species in the longwave region. More...
 
pure complex(kind=musica_dk) function, dimension(number_of_bands) net_shortwave_refractive_index (this, mode_state, number_of_bands)
 Calculates the net refractive index for the mode based on weighted contributions from constituent species in the shortwave region. More...
 
class(optics_t) function, pointer new_optics (this, property, output_grid, interpolation_strategy)
 Creates an optics_t object for a given property. More...
 
class(aerosol_state_t) function, pointer new_state (this)
 Creates a new state object for the mode. More...
 
real(kind=musica_dk) elemental function number_mixing_ratio__num_mol (this, mode_state)
 Returns the particle number mixing ratio per mole air [# mol-1]. More...
 
subroutine print_state (this, aerosol_state, io_unit)
 Outputs the current mode state. More...
 
subroutine randomize (this)
 Set the mode state to a random, but reasonable, state. For testing only. More...
 
integer function raw_size (this)
 Returns the number of doubles needed to hold the raw mode state. More...
 
type(wavelength_grid_t) function shortwave_grid (this)
 Returns the native shortwave optics wavelength grid. More...
 
subroutine shortwave_optics_array (this, environmental_state, aerosol_state, optics)
 Returns shortwave optical properties. More...
 
subroutine shortwave_optics_scalar (this, environmental_state, aerosol_state, optics)
 Returns shortwave optical properties. More...
 
pure real(kind=musica_dk) function, dimension(number_of_bands) specific_absorption__m2_kg (this, mode_state, number_of_bands, number_of_coefficients, coefficients, size_function, max_absorption)
 Calculates the specific absorption [m2 kg-1] for a given Chebyshev function. More...
 
pure real(kind=musica_dk) function, dimension(number_of_bands) specific_extinction__m2_kg (this, mode_state, number_of_bands, number_of_coefficients, coefficients, size_function, optics_lookup)
 Calculates the specific extinction [m2 kg-1] for a given Chebyshev function. More...
 
real(kind=musica_dk) elemental function wet_number_mode_diameter__m (this, mode_state)
 Returns the mode diameter [m] of the number distribution for the mode. More...
 
real(kind=musica_dk) elemental function wet_number_mode_radius__m (this, mode_state)
 Returns the mode radius [m] of the number distribution for the mode. More...
 
real(kind=musica_dk) elemental function wet_surface_mode_diameter__m (this, mode_state)
 Returns the mode diameter [m] of the surface area distribution for the mode. More...
 
real(kind=musica_dk) elemental function wet_surface_mode_radius__m (this, mode_state)
 Returns the mode radius [m] of the surface area distribution for the mode. More...
 
elemental real(kind=musica_dk) function wet_volume_to_mass_mixing_ratio__m3_kg (this, mode_state)
 Returns the wet volume to mass mixing ratio for the mode [m3 kg_dryair-1]. More...
 

Detailed Description

The mode_t type and related functions.

Function/Subroutine Documentation

◆ add_longwave_optics()

subroutine mam_mode::add_longwave_optics ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(mode_state_t), intent(in)  mode_state,
class(optics_ptr), dimension(:), intent(inout)  optics 
)

Calculates longwave optical properties for a given state and adds them to the set of optical properties.

Definition at line 434 of file mode.F90.

◆ add_shortwave_optics()

subroutine mam_mode::add_shortwave_optics ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(mode_state_t), intent(in)  mode_state,
class(optics_ptr), dimension(:), intent(inout)  optics 
)

Calculates shortwave optical properties for a given state and adds them to the set of optical properties.

Todo:
Are single-scatter albedo and layer_aod the correct terms to use?

Definition at line 353 of file mode.F90.

◆ asymmetry_factor()

pure real(kind=musica_dk) function, dimension( number_of_bands ) mam_mode::asymmetry_factor ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state,
integer, intent(in)  number_of_bands,
integer, intent(in)  number_of_coefficients,
real(kind=musica_dk), dimension( number_of_coefficients, number_of_bands ), intent(in)  coefficients,
real(kind=musica_dk), dimension( number_of_coefficients ), intent(in)  size_function 
)

Calculates the asymmetry parameter [unitless ratio] for a given Chebyshev function.

Returns
Asymmetry factor by wavelength
Parameters
[in]thisMAM mode
[in]mode_stateMAM mode state
[in]number_of_bandsNumber of wavelength bands
[in]number_of_coefficientsNumber of Chebyshev coefficients
[in]coefficientsChebyshev coefficients for asymmetry parameter calculation
[in]size_functionChebyshev function for the current surface mode radius

Definition at line 787 of file mode.F90.

◆ constructor()

type(mode_t) function mam_mode::constructor ( class(config_t), intent(inout)  config)
private

Constructor of mode_t objects.

Definition at line 87 of file mode.F90.

Referenced by mam_mode::mode_t::asymmetry_factor().

◆ dump_state()

subroutine mam_mode::dump_state ( class(mode_state_t), intent(in)  this,
real(kind=musica_dk), dimension(:), intent(inout)  raw_state,
integer, intent(inout), optional  index 
)
private

Dumps the aerosol state into the raw state array.

Definition at line 861 of file mode.F90.

◆ geometric_mean_diameter_of_number_distribution__m()

elemental real(kind=musica_dk) function mam_mode::geometric_mean_diameter_of_number_distribution__m ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)

Returns the geometric mean diameter [m] of the number distribution of the mode.

Definition at line 491 of file mode.F90.

◆ geometric_standard_deviation_of_number_distribution()

elemental real(kind=musica_dk) function mam_mode::geometric_standard_deviation_of_number_distribution ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the geometric standard deviation of the number distribution of the mode.

Definition at line 506 of file mode.F90.

◆ load_state()

subroutine mam_mode::load_state ( class(mode_state_t), intent(inout)  this,
real(kind=musica_dk), dimension(:), intent(in)  raw_state,
integer, intent(inout), optional  index 
)
private

Loads raw mode state data to the mode_state_t object.

Definition at line 839 of file mode.F90.

◆ longwave_grid()

type(wavelength_grid_t) function mam_mode::longwave_grid ( class(mode_t), intent(in)  this)

Returns the native longwave optics wavelength grid.

Definition at line 289 of file mode.F90.

Referenced by longwave_grid().

◆ longwave_optics_array()

subroutine mam_mode::longwave_optics_array ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(aerosol_state_t), intent(in)  aerosol_state,
type(optics_ptr), dimension(:), intent(inout)  optics 
)

Returns longwave optical properties.

Definition at line 244 of file mode.F90.

◆ longwave_optics_scalar()

subroutine mam_mode::longwave_optics_scalar ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(aerosol_state_t), intent(in)  aerosol_state,
class(optics_t), intent(inout), target  optics 
)

Returns longwave optical properties.

Definition at line 221 of file mode.F90.

◆ net_longwave_refractive_index()

pure complex(kind=musica_dk) function, dimension( number_of_bands ) mam_mode::net_longwave_refractive_index ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state,
integer, intent(in)  number_of_bands 
)
private

Calculates the net refractive index for the mode based on weighted contributions from constituent species in the longwave region.

Definition at line 650 of file mode.F90.

◆ net_shortwave_refractive_index()

pure complex(kind=musica_dk) function, dimension( number_of_bands ) mam_mode::net_shortwave_refractive_index ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state,
integer, intent(in)  number_of_bands 
)
private

Calculates the net refractive index for the mode based on weighted contributions from constituent species in the shortwave region.

Definition at line 616 of file mode.F90.

◆ new_optics()

class(optics_t) function, pointer mam_mode::new_optics ( class(mode_t), intent(in)  this,
class(property_t), intent(in)  property,
type(wavelength_grid_t), intent(in)  output_grid,
procedure(interpolation_strategy_i), optional  interpolation_strategy 
)

Creates an optics_t object for a given property.

Definition at line 143 of file mode.F90.

Referenced by new_optics().

◆ new_state()

class(aerosol_state_t) function, pointer mam_mode::new_state ( class(mode_t), intent(in)  this)

Creates a new state object for the mode.

Definition at line 125 of file mode.F90.

Referenced by new_state().

◆ number_mixing_ratio__num_mol()

real(kind=musica_dk) elemental function mam_mode::number_mixing_ratio__num_mol ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the particle number mixing ratio per mole air [# mol-1].

Todo:
is this per mole dry or wet air?

Definition at line 602 of file mode.F90.

◆ print_state()

subroutine mam_mode::print_state ( class(mode_t), intent(in)  this,
class(aerosol_state_t), intent(in)  aerosol_state,
integer, intent(in), optional  io_unit 
)

Outputs the current mode state.

Parameters
[in]thisMAM mode
[in]aerosol_stateMAM mode state
[in]io_unitOptional output unit (defaults to 6)

Definition at line 303 of file mode.F90.

◆ randomize()

subroutine mam_mode::randomize ( class(mode_state_t), intent(inout)  this)
private

Set the mode state to a random, but reasonable, state. For testing only.

Definition at line 883 of file mode.F90.

◆ raw_size()

integer function mam_mode::raw_size ( class(mode_state_t), intent(in)  this)

Returns the number of doubles needed to hold the raw mode state.

Definition at line 828 of file mode.F90.

Referenced by raw_size().

◆ shortwave_grid()

type(wavelength_grid_t) function mam_mode::shortwave_grid ( class(mode_t), intent(in)  this)

Returns the native shortwave optics wavelength grid.

Definition at line 275 of file mode.F90.

Referenced by shortwave_grid().

◆ shortwave_optics_array()

subroutine mam_mode::shortwave_optics_array ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(aerosol_state_t), intent(in)  aerosol_state,
type(optics_ptr), dimension(:), intent(inout)  optics 
)

Returns shortwave optical properties.

Definition at line 190 of file mode.F90.

◆ shortwave_optics_scalar()

subroutine mam_mode::shortwave_optics_scalar ( class(mode_t), intent(in)  this,
class(environmental_state_t), intent(in)  environmental_state,
class(aerosol_state_t), intent(in)  aerosol_state,
class(optics_t), intent(inout), target  optics 
)

Returns shortwave optical properties.

Definition at line 167 of file mode.F90.

◆ specific_absorption__m2_kg()

pure real(kind=musica_dk) function, dimension( number_of_bands ) mam_mode::specific_absorption__m2_kg ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state,
integer, intent(in)  number_of_bands,
integer, intent(in)  number_of_coefficients,
real(kind=musica_dk), dimension( number_of_coefficients, number_of_bands ), intent(in)  coefficients,
real(kind=musica_dk), dimension( number_of_coefficients ), intent(in)  size_function,
real(kind=musica_dk), dimension( number_of_bands ), intent(in), optional  max_absorption 
)
private

Calculates the specific absorption [m2 kg-1] for a given Chebyshev function.

Todo:
is "specific absorption" (m2 kg-1) the correct name/units for what is returned from this function?
Parameters
[in]coefficientsChebyshev coefficients for absorption calculation
[in]size_functionChebyshev function for the current surface mode radius
[in]max_absorptionMaximum values for specific absorption for each band [m2 kg-1]

Definition at line 687 of file mode.F90.

Referenced by specific_absorption__m2_kg().

◆ specific_extinction__m2_kg()

pure real(kind=musica_dk) function, dimension( number_of_bands ) mam_mode::specific_extinction__m2_kg ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state,
integer, intent(in)  number_of_bands,
integer, intent(in)  number_of_coefficients,
real(kind=musica_dk), dimension( number_of_coefficients, number_of_bands ), intent(in)  coefficients,
real(kind=musica_dk), dimension( number_of_coefficients ), intent(in)  size_function,
class(optics_lookup_t), intent(in)  optics_lookup 
)

Calculates the specific extinction [m2 kg-1] for a given Chebyshev function.

Todo:
is "specific extinction" (m2 kg-1) the correct name/units for what is returned from this function?
Returns
Specific extinction by wavelength
Parameters
[in]thisMAM mode
[in]mode_stateMAM mode state
[in]number_of_bandsNumber of wavelength bands
[in]number_of_coefficientsNumber of Chebyshev coefficients
[in]coefficientsChebyshev coefficients for extinction calculation
[in]size_functionChebyshev function for the current surface mode radius
[in]optics_lookupOptics lookup table
Todo:
is the data returned from the lookup table in ln( m2/kg )?
Todo:
in the original code, this is labelled as converting "m2/kg water to m2/kg aerosol", but it appears to convert instead to m2/kg dry air - is this right???

Definition at line 732 of file mode.F90.

Referenced by specific_extinction__m2_kg().

◆ wet_number_mode_diameter__m()

real(kind=musica_dk) elemental function mam_mode::wet_number_mode_diameter__m ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the mode diameter [m] of the number distribution for the mode.

Definition at line 556 of file mode.F90.

◆ wet_number_mode_radius__m()

real(kind=musica_dk) elemental function mam_mode::wet_number_mode_radius__m ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the mode radius [m] of the number distribution for the mode.

Definition at line 543 of file mode.F90.

◆ wet_surface_mode_diameter__m()

real(kind=musica_dk) elemental function mam_mode::wet_surface_mode_diameter__m ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the mode diameter [m] of the surface area distribution for the mode.

Definition at line 570 of file mode.F90.

◆ wet_surface_mode_radius__m()

real(kind=musica_dk) elemental function mam_mode::wet_surface_mode_radius__m ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the mode radius [m] of the surface area distribution for the mode.

Definition at line 583 of file mode.F90.

◆ wet_volume_to_mass_mixing_ratio__m3_kg()

elemental real(kind=musica_dk) function mam_mode::wet_volume_to_mass_mixing_ratio__m3_kg ( class(mode_t), intent(in)  this,
class(mode_state_t), intent(in)  mode_state 
)
private

Returns the wet volume to mass mixing ratio for the mode [m3 kg_dryair-1].

Definition at line 520 of file mode.F90.