This module is used for outputing model results in a list of stations (not gridded arrays). The user needs to supply a list of stations with lat, lon values of each station. Data at a single point (I,J) that is closest to each station will be written to a file. No interpolation is made when a station is between two or more grid points.
In the output file, a 3D field will have a format of array(n1,n2) and a 2D field is array(n1) where n1 is number of stations and n2 is number of vertical levels or depths.
More...
|
subroutine, public | station_data_init () |
| Read in lat and lon of each station
Create station_id based on lat, lon
Read station_data_table, initialize output_fields and output files. More...
|
|
subroutine | check_duplicate_output_fields () |
| check_duplicate_output_fields takes the data pairs (output_name and output_file) and (module_name and output_name) and ensures each pair is unique. More...
|
|
integer function | get_station_id (lat, lon) |
| get_station_id is passed the station's distinct lat and lon to determine what the station's id is. More...
|
|
subroutine | init_file (filename, output_freq, output_units, time_units) |
| init_file initializes files to write station data to. Data is formatted as number of units of measurement at X output frequency over amount of time dt More...
|
|
subroutine | init_output_field (module_name, field_name, file_name, time_method, pack) |
| init_output_field initializes output_field attributes More...
|
|
integer function | find_file (name) |
| find_file finds the index of a requested file More...
|
|
integer function | register_station_field2d (module_name, fieldname, glo_lat, glo_lon, init_time, domain, longname, units) |
| register_station_field2d registers a new station field with user input.
register_station_field2d is overloaded by register_station_field3d as a part of the register_station_field interface. More...
|
|
integer function | register_station_field3d (module_name, fieldname, glo_lat, glo_lon, levels, init_time, domain, longname, units) |
| register_station_field3d registers a new station field with user input.
register_station_field3d is overloaded by register_station_field2d as a part of the register_station_field interface. More...
|
|
integer function | find_output_field (module_name, field_name) |
| find_output_field returns the index of the reuqested field within the output_fields array. More...
|
|
subroutine | opening_file (file) |
| opening_file opens a file and writes axis meta_data for all files (only on ROOT PE, do nothing on other PEs) More...
|
|
subroutine | send_station_data_2d (field_id, data, time) |
| send_station_data_2d sends data to the root PE, which then sends data to staton_data_out to be sent to files.
send_station_data_2d is overloaded by send_station_data_3d as a part of the send_station_data interface. More...
|
|
subroutine | send_station_data_3d (field_id, data, time) |
| send_station_data_3d sends data to the root PE, which then sends data to staton_data_out to be sent to files.
send_station_data_3d is overloaded by send_station_data_2d as a part of the send_station_data interface. More...
|
|
subroutine | station_data_out (file, field, data, time, final_call_in) |
| station_data_out is responsible for sending station data to files. More...
|
|
subroutine, public | station_data_end (time) |
| Must be called after the last time step to write the buffer content. More...
|
|
|
integer, parameter | max_fields_per_file = 150 |
|
integer, parameter | max_files = 31 |
|
integer | num_files = 0 |
|
integer | num_stations = 0 |
|
integer | max_stations = 20 |
|
integer | max_output_fields = 100 |
|
integer | num_output_fields = 0 |
|
real | empty = 0.0 |
|
real | missing = 1.E20 |
|
logical | module_is_initialized = .false. |
|
logical | need_write_axis = .true. |
|
integer | base_year |
|
integer | base_month |
|
integer | base_day |
|
integer | base_hour |
|
integer | base_minute |
|
integer | base_second |
|
type(time_type) | base_time |
|
character(len=10), dimension(6) | time_unit_list = (/'seconds ', 'minutes ', 'hours ', 'days ', 'months ', 'years '/) |
|
integer, parameter | every_time = 0 |
|
integer, parameter | end_of_run = -1 |
|
character(len=256) | global_descriptor |
|
character(len=7) | avg_name = 'average' |
|
integer | total_pe |
|
integer | lat_axis |
|
integer | lon_axis |
|
integer, dimension(:), allocatable | pelist |
|
character(len=32) | pelist_name |
|
type(global_field_type), save | global_field |
|
type(file_type), dimension(max_files), save | files |
|
type(group_field_type), dimension(:), allocatable, save | output_fields |
|
type(station_type), dimension(:), allocatable | stations |
|
type(diag_fieldtype), save | diag_field |
|
This module is used for outputing model results in a list of stations (not gridded arrays). The user needs to supply a list of stations with lat, lon values of each station. Data at a single point (I,J) that is closest to each station will be written to a file. No interpolation is made when a station is between two or more grid points.
In the output file, a 3D field will have a format of array(n1,n2) and a 2D field is array(n1) where n1 is number of stations and n2 is number of vertical levels or depths.
- Author
- Giang Nong Giang.nosp@m..Non.nosp@m.g@gfd.nosp@m.l.no.nosp@m.aa.go.nosp@m.v
Here are some basic steps of how to use station_data_mod:
-
Call
data_station_init
User needs to supply 2 tables: list_stations
and station_data_table
as follows:
Example of list_stations : | Example of station_data_table : |
station_id | lat | lon |
station_1 | 20.4 | 100.8 |
| | |
| | |
| | |
| | |
| | |
| | |
| | |
| | |
| | |
|
General descriptor: |
Am2p14 station data |
Start time (should be the same as model's initial time): |
19800101 |
File information: |
filename | output_frequency | frequency_unit | time_axis_unit |
"ocean_day" | 1 | "days" | "hours" |
|
Field information: |
module | field_name | filename | time_method | pack |
Ice_mod | temperature | ocean_day | .true. | 2 |
Ice_mod | pressure | ocean_day | .false. | 2 |
|
|
-
Call
register_station_field
to register each field that needs to be written to a file, the call register_station_field
returns a field_id that will be used later in send_station_data
-
Call
send_station_data
will send data at each station in the list to a file
-
Finally, call
station_data_end
after the last time step.
Modules Included:
Module Name | Functions Included |
axis_utils_mod | nearest_index |
mpp_io_mod | mpp_open, MPP_RDONLY, MPP_ASCII, mpp_close, MPP_OVERRWR, MPP_NETCDF, mpp_write_meta, MPP_SINGLE, mpp_write, fieldtype, mpp_flush |
fms_mod | error_mesg, FATAL, WARNING, stdlog, write_version_number, mpp_pe, lowercase, stdout, close_file, open_namelist_file, check_nml_error |
mpp_mod | mpp_npes, mpp_sync, mpp_root_pe, mpp_send, mpp_recv, mpp_max, mpp_get_current_pelist, input_nml_file, COMM_TAG_1, COMM_TAG_2, COMM_TAG_3, COMM_TAG_4 |
mpp_domains_mod | domain2d, mpp_get_compute_domain |
diag_axis_mod | diag_axis_init |
diag_output_mod | write_axis_meta_data, write_field_meta_data, diag_fieldtype, done_meta_data |
diag_manager_mod | get_date_dif, DIAG_SECONDS, DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS, DIAG_MONTHS, DIAG_YEARS |
diag_util_mod | diag_time_inc |
time_manager_mod | operator(>), operator(>=), time_type, get_calendar_type, NO_CALENDAR, set_time set_date, increment_date, increment_time |
integer function station_data_mod::register_station_field3d |
( |
character(len=*), intent(in) |
module_name, |
|
|
character(len=*), intent(in) |
fieldname, |
|
|
real, dimension(:), intent(in) |
glo_lat, |
|
|
real, dimension(:), intent(in) |
glo_lon, |
|
|
real, dimension(:), intent(in) |
levels, |
|
|
type(time_type), intent(in) |
init_time, |
|
|
type(domain2d), intent(in) |
domain, |
|
|
character(len=*), intent(in), optional |
longname, |
|
|
character(len=*), intent(in), optional |
units |
|
) |
| |
|
private |
register_station_field3d registers a new station field with user input.
register_station_field3d is overloaded by register_station_field2d as a part of the register_station_field interface.
- Exceptions
-
FATAL,outside global latitude values | |
FATAL,outside global longitude values | |
FATAL,Error in global index of station | |
WARNING,<module_name>/<field_name> NOT found in station_data table | |
FATAL,Error in determining local_num_station | |
- Returns
- write field meta data on ROOT PE only allocate buffer
- Parameters
-
[in] | module_name | Module name of the station being registered |
[in] | fieldname | Field name of the tation being registered |
[in] | glo_lat | Global latitude location of the station (in X direction) |
[in] | glo_lon | Global longitude location of the station (in Y direction) |
[in] | levels | Global elevation location of the station (in Z direction) |
[in] | domain | Optional 2d domain representation |
[in] | init_time | Initial time of the station |
[in] | longname | Optional longname to identify station |
[in] | units | Optional units for station measurement data |
determine global index of this field in all stations
determine local index of this field in all stations , local index starts from 1
get the position of this field in the array output_fields
fill out list of available stations in this PE
deal with axes
Definition at line 761 of file station_data.F90.
subroutine station_data_mod::send_station_data_3d |
( |
integer, intent(in) |
field_id, |
|
|
real, dimension(:,:,:), intent(in) |
data, |
|
|
type(time_type), intent(in) |
time |
|
) |
| |
|
private |
send_station_data_3d sends data to the root PE, which then sends data to staton_data_out to be sent to files.
send_station_data_3d is overloaded by send_station_data_2d as a part of the send_station_data interface.
- Exceptions
-
FATAL,Station data NOT initialized | |
FATAL,counter=0 for averaged field | |
FATAL,Global field contains MISSING field | |
FATAL,Local index out of range for field | |
- Parameters
-
[in] | field_id | Field id of the station data being recorded |
[in] | data | Data of the station being recorded |
[in] | time | Time of the station being recorded |
compare time with next_output
get max_counter if the field is averaged
receive local data from all PEs
send global_buffer content to file
check if global_field contains any missing values
clear buffer, increment next_output time and reset counter on ALL PEs
accumulate buffer only
Definition at line 1020 of file station_data.F90.
subroutine, public station_data_mod::station_data_end |
( |
type(time_type), intent(in) |
time | ) |
|
Must be called after the last time step
to write the buffer content.
- Exceptions
-
FATAL,counter=0 for averaged field | |
FATAL,Global field contains MISSING field | |
- Parameters
-
ALL PEs, including root PE, must send data to root PE
get max_counter if the field is averaged
only root PE receives local data from all PEs
send global_buffer content to file
check if global_field contains any missing values
deallocate field buffer
Definition at line 1171 of file station_data.F90.
subroutine, public station_data_mod::station_data_init |
( |
| ) |
|
Read in lat and lon of each station
Create station_id based on lat, lon
Read station_data_table, initialize output_fields and output files.
read namelist
read list of stations
read station_data table
Read in the global file labeling string
Read in the base date
If no calendar, ignore year and month
Definition at line 395 of file station_data.F90.