MPAS-JEDI
|
Data Types | |
type | mpasjedi_getvalues_base |
type | mpasjedi_getvalues |
Functions/Subroutines | |
subroutine, public | getvalues_base_create (self, geom, locs, f_conf) |
Linked list implementation. More... | |
subroutine | create (self, geom, locs, f_conf) |
GetValues class 'create' logic. More... | |
subroutine, public | getvalues_base_delete (self) |
GetValues base class 'delete' logic. More... | |
subroutine | delete (self) |
GetValues class 'delete' logic. More... | |
subroutine, public | fill_geovals (self, geom, state, t1, t2, locs, gom) |
Interpolates from geovar mpas_fields to populate ufo_geovals. More... | |
subroutine | initialize_uns_interp (self, grid, lats_obs, lons_obs) |
Initializes an unstructured interpolation type. More... | |
subroutine | integer_interpolation_bump (self, ngrid, nlocs, data_in, obs_field_int, gom, jvar, time_mask) |
Performs interpolation of integer fields using BUMP. More... | |
subroutine | integer_interpolation_unstructured (self, ngrid, nlocs, data_in, work_field, gom, jvar, time_mask) |
Performs interpolation of integer fields using unstructured interpolation. More... | |
Variables | |
type(registry_t), public | mpas_getvalues_registry |
Linked list interface - defines registry_t type. More... | |
character(len=1024) | message |
|
private |
GetValues class 'create' logic.
create This subroutine contstructs an mpasjedi_getvalues object class instance.
[in,out] | self | getvalues self |
[in] | geom | geometry (mpas mesh) |
[in] | locs | ufo geovals (obs) locations |
[in] | f_conf | configuration |
Definition at line 147 of file mpasjedi_getvalues_mod.F90.
|
private |
GetValues class 'delete' logic.
delete This subroutine deletes (frees memory) for an mpasjedi_getvalues object class instance.
[in,out] | self | getvalues self |
Definition at line 178 of file mpasjedi_getvalues_mod.F90.
subroutine, public mpasjedi_getvalues_mod::fill_geovals | ( | class(mpasjedi_getvalues_base), intent(inout) | self, |
type(mpas_geom), intent(in) | geom, | ||
type(mpas_fields), intent(in) | state, | ||
type(datetime), intent(in) | t1, | ||
type(datetime), intent(in) | t2, | ||
type(ufo_locations), intent(in) | locs, | ||
type(ufo_geovals), intent(inout) | gom | ||
) |
Interpolates from geovar mpas_fields to populate ufo_geovals.
fill_geovals This subroutine populates the variables in a ufo_geovals object by interpolating the state variables in an mpas_fields object. This is the non-linear subroutine used in both GetValues and LinearGetValues classes
[in,out] | self | getvalues_base self |
[in] | geom | geometry (mpas mesh) |
[in] | state | state containing geovars |
[in] | t1 | time window begin |
[in] | t2 | time window end |
[in] | locs | observation locations |
[in,out] | gom | geovals |
Definition at line 190 of file mpasjedi_getvalues_mod.F90.
subroutine, public mpasjedi_getvalues_mod::getvalues_base_create | ( | class(mpasjedi_getvalues_base), intent(inout) | self, |
type(mpas_geom), intent(in) | geom, | ||
type(ufo_locations), intent(in) | locs, | ||
type(fckit_configuration), intent(in) | f_conf | ||
) |
Linked list implementation.
GetValues base class 'create' logic
getvalues_base_create This subroutine populates the getvalues_base class members. This subroutine is called from the 'create' subroutines of all derived classes. (i.e. getvalues and lineargetvalues)
[in,out] | self | getvalues_base self |
[in] | geom | geometry (mpas mesh) |
[in] | locs | ufo geovals (obs) locations |
[in] | f_conf | configuration |
Definition at line 99 of file mpasjedi_getvalues_mod.F90.
subroutine, public mpasjedi_getvalues_mod::getvalues_base_delete | ( | class(mpasjedi_getvalues_base), intent(inout) | self | ) |
GetValues base class 'delete' logic.
getvalues_base_delete This subroutine deletes (frees memory) for the getvalues_base class members. This subroutine is called from the 'delete' subroutines of all derived classes. (i.e. getvalues and lineargetvalues)
[in,out] | self | getvalues_base self |
Definition at line 163 of file mpasjedi_getvalues_mod.F90.
|
private |
Initializes an unstructured interpolation type.
initialize_uns_interp This subroutine calls unsinterpcreate, which calculates the barycentric weights used to interpolate data between the mpas mesh locations (grid) and the observation locations.
[in] | grid | mpas mesh data |
[in] | lats_obs | latitudes of obs |
[in] | lons_obs | longitudes of obs |
Definition at line 319 of file mpasjedi_getvalues_mod.F90.
|
private |
Performs interpolation of integer fields using BUMP.
integer_interpolation_bump This subroutine performs the interpolation of integer-valued fields (i.e. types) using BUMP
[in] | ngrid | number of cells in model mesh |
[in] | nlocs | number of locations for obs |
[in] | data_in | data to interpolate |
[in,out] | obs_field_int | output array of interpolated data |
[in,out] | gom | output geoVaLs |
[in] | jvar | gom % geovals index |
[in] | time_mask | mask for time window |
Definition at line 362 of file mpasjedi_getvalues_mod.F90.
|
private |
Performs interpolation of integer fields using unstructured interpolation.
integer_interpolation_unstructured This subroutine performs the interpolation of integer-valued fields (i.e. types) using unstructured interpolation
[in] | ngrid | number of cells in model mesh |
[in] | nlocs | number of locations for obs |
[in] | data_in | data to interpolate |
[in,out] | work_field | (ngrid,1) array |
[in,out] | gom | output geoVaLs |
[in] | jvar | gom % geovals index |
[in] | time_mask | mask for time window |
Definition at line 397 of file mpasjedi_getvalues_mod.F90.
|
private |
Definition at line 83 of file mpasjedi_getvalues_mod.F90.
type(registry_t), public mpasjedi_getvalues_mod::mpas_getvalues_registry |
Linked list interface - defines registry_t type.
Global registry
Definition at line 79 of file mpasjedi_getvalues_mod.F90.