mam  v1.0
A Modal Aerosol Model
species.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_species module
6
7!> The species_t type and related functions
9
10 use musica_constants, only : musica_dk
11 use musica_string, only : string_t
12
13 implicit none
14 private
15
16 public :: species_t
17
18 !> \todo species should read in the wavelength grid for refractive
19 !! indices instead of assuming it matches the MAM grid
20 type :: species_t
21 private
22 type(string_t) :: name_
23 real(kind=musica_dk) :: density__kg_m3_
24 complex(kind=musica_dk), allocatable :: shortwave_refractive_index_(:)
25 complex(kind=musica_dk), allocatable :: longwave_refractive_index_(:)
26 contains
27 procedure :: name
28 procedure :: volume__m3
31 end type species_t
32
33 interface species_t
34 module procedure :: constructor
35 end interface
36
37contains
38
39!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
40
41 !> Constructs species_t objects
42 !!
43 !! \todo The complex refractive indices in the NetCDF files use a different
44 !! convention for specifying the complex index of refraction that
45 !! results in the need to take the absolute value of the imaginary
46 !! part to get the complex index of refraction used in MAM.
47 !! The NetCDF files should be updated to use positive values for the
48 !! imaginary part.
49 function constructor( config ) result( new_obj )
50
51 use musica_assert, only : assert_msg
52 use musica_config, only : config_t
53 use musica_io, only : io_t
54 use musica_io_netcdf, only : io_netcdf_t
55
56 type(species_t) :: new_obj
57 class(config_t), intent(inout) :: config
58
59 character(len=*), parameter :: my_name = "species_t constructor"
60 class(io_t), pointer :: file_ptr
61 type(string_t) :: file_path
62 type(string_t) :: sw_real_ri, sw_imag_ri, lw_real_ri, lw_imag_ri
63 type(config_t) :: optics
64 real(kind=musica_dk), allocatable :: real_values(:), imag_values(:)
65 logical :: found
66
67 call config%get( "name", new_obj%name_, my_name )
68 call config%get( "density", "kg m-3", new_obj%density__kg_m3_, my_name )
69 call config%get( "optics", optics, my_name, found = found )
70 if( found ) then
71 call optics%get( "file path", file_path, my_name )
72 file_ptr => io_netcdf_t( file_path )
73 call optics%get( "real refractive index - shortwave", sw_real_ri, &
74 my_name )
75 call optics%get( "real refractive index - longwave", lw_real_ri, &
76 my_name )
77 call optics%get( "imaginary refractive index - shortwave", sw_imag_ri, &
78 my_name )
79 call optics%get( "imaginary refractive index - longwave", lw_imag_ri, &
80 my_name )
81 call file_ptr%read( sw_real_ri, real_values, my_name )
82 call file_ptr%read( sw_imag_ri, imag_values, my_name )
83 call assert_msg( 434735228, &
84 size( real_values ) .eq. size( imag_values ), &
85 "Length mismatch for shortwave refractive index "// &
86 "arrays for species '"//new_obj%name_//"'" )
87 allocate( new_obj%shortwave_refractive_index_( size( real_values ) ) )
88 new_obj%shortwave_refractive_index_(:) = &
89 cmplx( real_values(:), abs( imag_values(:) ), kind = musica_dk )
90 deallocate( real_values )
91 deallocate( imag_values )
92 call file_ptr%read( lw_real_ri, real_values, my_name )
93 call file_ptr%read( lw_imag_ri, imag_values, my_name )
94 call assert_msg( 979080793, &
95 size( real_values ) .eq. size( imag_values ), &
96 "Length mismatch for longwave refractive index "// &
97 "arrays for species '"//new_obj%name_//"'" )
98 allocate( new_obj%longwave_refractive_index_( size( real_values ) ) )
99 new_obj%longwave_refractive_index_(:) = &
100 cmplx( real_values(:), abs( imag_values(:) ), kind = musica_dk )
101 deallocate( file_ptr )
102 else
103 new_obj%shortwave_refractive_index_(:) = cmplx( 0.0, 0.0 )
104 new_obj%longwave_refractive_index_(:) = cmplx( 0.0, 0.0 )
105 end if
106
107 end function constructor
108
109!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
110
111 !> Returns the species name
112 type(string_t) function name( this )
113
114 class(species_t), intent(in) :: this
115
116 name = this%name_
117
118 end function name
119
120!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
121
122 !> Returns the total volume occupied by the species for a given mass
123 real(kind=musica_dk) elemental function volume__m3( this, species_mass__kg )
124
125 class(species_t), intent(in) :: this
126 real(kind=musica_dk), intent(in) :: species_mass__kg
127
128 volume__m3 = species_mass__kg / this%density__kg_m3_
129
130 end function volume__m3
131
132!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
133
134 !> Returns the shortwave refractive index at a specified band index
135 complex(kind=musica_dk) elemental function shortwave_refractive_index( this,&
136 band )
137
138 class(species_t), intent(in) :: this
139 integer, intent(in) :: band
140
141 shortwave_refractive_index = this%shortwave_refractive_index_( band )
142
143 end function shortwave_refractive_index
144
145!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
146
147 !> Returns the longwave refractive index at a specified band index
148 complex(kind=musica_dk) elemental function longwave_refractive_index( this, &
149 band )
150
151 class(species_t), intent(in) :: this
152 integer, intent(in) :: band
153
154 longwave_refractive_index = this%longwave_refractive_index_( band )
155
156 end function longwave_refractive_index
157
158!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
159
160end module mam_species
The species_t type and related functions.
Definition: species.F90:8
type(species_t) function constructor(config)
Constructs species_t objects.
Definition: species.F90:50
type(string_t) function name(this)
Returns the species name.
Definition: species.F90:113
real(kind=musica_dk) elemental function volume__m3(this, species_mass__kg)
Returns the total volume occupied by the species for a given mass.
Definition: species.F90:124
complex(kind=musica_dk) elemental function longwave_refractive_index(this, band)
Returns the longwave refractive index at a specified band index.
Definition: species.F90:150
complex(kind=musica_dk) elemental function shortwave_refractive_index(this, band)
Returns the shortwave refractive index at a specified band index.
Definition: species.F90:137