obsSubSpaceData_mod

link to source code

Dependency Diagrams:

obsSubSpaceData_mod.svg

Direct Dependency Diagram

obsSubSpaceData_mod_rev.svg

Reverse Dependency Diagram

Description

MODULE obsSubSpaceData_mod (prefix=’oss’ category=’6. High-level data objects’)

Purpose

Repository of obs space structures, arrays, and routines specific to obs data pertinent to subspaces of the overall ObsSpaceData. Most tools are independent of ObsSpaceData and can be used by themselves for the users’ application(s).

Quick access

Variables

oss_obsdata_get_header_code, struct_oss_obsdata

Routines

obsdata_extra_code_test(), obsdata_get_header_code_i(), obsdata_get_header_code_r(), obsdata_set_index(), oss_comboidlist(), oss_get_comboidlist(), oss_obsdata_add_data1d(), oss_obsdata_alloc(), oss_obsdata_code_len(), oss_obsdata_dealloc(), oss_obsdata_get_array1d(), oss_obsdata_get_array2d(), oss_obsdata_get_data1d(), oss_obsdata_get_element(), oss_obsdata_mpiallgather()

Needed modules

  • codeprecision_mod: MODULE codePrecision_mod (prefix=’pre’ category=’8. Low-level utilities and constants’)

  • utilities_mod: MODULE utilities_mod (prefix=’utl’ category=’8. Low-level utilities and constants’)

  • mathphysconstants_mod: MODULE mathPhysConstants_mod (prefix=’mpc’ category=’8. Low-level utilities and constants’)

  • midasmpi_mod: MODULE midasMpi_mod (prefix=’mmpi’ category=’8. Low-level utilities and constants’)

  • bufr_mod: MODULE bufr_mod (prefix=’bufr’ category=’8. Low-level utilities and constants’)

  • obsspacedata_mod: MODULE obsSpaceData_mod (prefix=’obs’ category=’6. High-level data objects’)

Types

  • type  obssubspacedata_mod/unknown_type
    Type fields
    • % code (*) [character ,pointer]

    • % data1d (*,*) [real ,pointer]

    • % data2d (*,*,*) [real ,pointer]

    • % dim1 [integer ]

    • % dim2 [integer ]

    • % irep [integer ]

    • % ndim [integer ]

    • % nrep [integer ]

Variables

  • obssubspacedata_mod/oss_obsdata_get_header_code [public]
  • obssubspacedata_mod/struct_oss_obsdata [public]

Subroutines and functions

subroutine  obssubspacedata_mod/oss_obsdata_alloc(obsdata, nrep, dim1[, dim2_opt])
Purpose

Allocates memory for structure struct_oss_obsdata to hold obs data file information. If dim2 is specified, then the data array associated with each observation/report will be 2D array. will be a 1D array if dim2 is not specified.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: data structure to allocate

  • nrep [integer ,in] :: max number of associated observations/reports for data

  • dim1 [integer ,in] :: first dimension length of the array associated to each observation/report :dim2: second dimension length of the array associated to each observation/report (optional)

Options

dim2_opt [integer ,in,]

Called from

bcsc_addbgstddev(), brpf_obssub_read(), obsf_obssub_read(), oopc_setupch(), oss_obsdata_add_data1d(), oss_obsdata_mpiallgather()

subroutine  obssubspacedata_mod/oss_obsdata_dealloc(obsdata)
Purpose

Deallocates memory for structure struct_oss_obsdata

Arguments

obsdata [struct_oss_obsdata ,inout]

Called from

chm_dealloc_obs_err_stddev(), oopc_deallocavgkern(), osd_readsqrthpht(), oss_obsdata_mpiallgather()

function  obssubspacedata_mod/oss_obsdata_get_element(obsdata, code, idim1[, stat_opt])
Purpose

Returns element of array from obsdata 1D data array. The returned element is the one with the specified identifying code.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: struct_oss_obsdata instance

  • code [character ,in] :: unique identifying code

  • idim1 [integer ,in] :: position of element in the dim1 axis

Return

element [real ] :: retrieved element from obsdata%data1d :stat: status code

Options

stat_opt [integer ,out,]

Called from

brpf_obssub_update(), chm_get_obs_err_stddev()

Call to

obsdata_set_index()

function  obssubspacedata_mod/oss_obsdata_get_array1d(obsdata, code[, stat_opt])
Purpose

Returns 1D data array from obsdata. The returned array is the one with the specified identifying code.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: struct_oss_obsdata instance

  • code [character ,in] :: unique identifying code

Return

array (obsdata%dim1) [real ] :: retrieved array from obsdata%data1d of dimension obsdata%dim1 :stat: search success (0 - found; 1 = no data; 2 = not found)

Options

stat_opt [integer ,out,]

Called from

brpf_obssub_update(), oss_obsdata_get_data1d()

Call to

obsdata_set_index()

function  obssubspacedata_mod/oss_obsdata_get_data1d(obsdata, lon, lat, date, time, stnid[, stat_opt])
Purpose

Extract 1D data array from structure according to input (lat,long,date,time,stnid)

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: Structure from which data elements are to be found

  • lon [real ,in] :: longitude real (radians)

  • lat [real ,in] :: latitude real (radians)

  • date [integer ,in] :: date in YYYYMMDD

  • time [integer ,in] :: time in HHMM

  • stnid [character ,in] :: station ID

Return

array (obsdata%dim1) [real ] :: Identified 1D array :stat: search success (0 - found; 1 = no data; 2 = not found)

Options

stat_opt [integer ,out,]

Called from

osd_readsqrthpht()

Call to

oss_obsdata_get_array1d()

function  obssubspacedata_mod/oss_obsdata_get_array2d(obsdata, code[, stat_opt])
Purpose

Returns 2D data array from obsdata. The returned array is the one with the specified identifying code.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: struct_oss_obsdata instance

  • code [character ,in] :: unique identifying code

Return

array (obsdata%dim1,obsdata%dim2) [real ] :: retrieved array from obsdata%data2d of dimension (obsdata%dim1,obsdata%dim2) :stat: search success (0 - found; 1 = no data; 2 = not found)

Options

stat_opt [integer ,out,]

Called from

bcsc_retrievebgstddev(), brpf_obssub_update(), oopc_getavgkern()

Call to

obsdata_set_index()

subroutine  obssubspacedata_mod/obsdata_set_index(obsdata, code[, stat_opt])
Purpose

Sets the position variable (irep) in struct_oss_obsdata to reference the record that matches the input identifying code.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: irep: current index position for the observation/report :index: obs index

  • code [character ,in] :: code for comparison to those in obsdata

  • obsdata :: irep: current index position for the observation/report

Options

stat_opt [integer ,out,] :: status of call (optional) :0: no errors :1: no reports available :2: report not found

Comments

If the optional argument stat_opt is provided and an error occurs, the error code will be returned and the abort will not be called to allow for error handling.

Called from

oss_obsdata_get_element(), oss_obsdata_get_array1d(), oss_obsdata_get_array2d()

Call to

utl_abort(), obsdata_extra_code_test()

function  obssubspacedata_mod/obsdata_extra_code_test(test_code, ref_code, ref_lat)
Purpose

Test matching of code values accounting for rare differences in stored lat (and lon) value(s) when codes are stored as strings in the form LAT–LON–YYYYMMDDHHMM* (ie. with >= oss_code_sublen characters).

Caveat

The current version assumes the only source of difference would stem from a shift to the nearest latitude of the analysis grid from points near the pole. (this source of difference identified by M. Sitwell) Also currently assumes that most poleward analysis grid latitudes are within 1 degree away from a pole.

Arguments

Arguments
  • test_code [character ,in] :: code for comparison to ref_code

  • ref_lat [integer ,in] :: latitude (x100) part of reference code

  • ref_code [character ,in] :: reference code

Return

found [logical ] :: logical indicating if a match has been found.

Called from

obsdata_set_index()

function  obssubspacedata_mod/obsdata_get_header_code_i(ilon, ilat, date, time, stnid)
Purpose

Generates a string code to identify an obervation by the header information in a BURP report. The BURP header information is saved as a string in the form LAT–LON–YYYYMMDDHHMMSTNID—-. Intention of this function is to be used for setting the unique identifier ‘code’ in struct_oss_obsdata. Can be called under the interface oss_obsdata_get_header_code.

Arguments

Arguments
  • ilon [integer ,in] :: longitude integer

  • ilat [integer ,in] :: latitude integer

  • date [integer ,in] :: date in YYYYMMDD

  • time [integer ,in] :: time in HHMM

  • stnid [character ,in] :: station ID

Return

code [character ] :: unique code

Called from

obsdata_get_header_code_r()

function  obssubspacedata_mod/obsdata_get_header_code_r(lon, lat, date, time, stnid)
Purpose

Generates a string code to identify an obervation by the header information in a BURP report. The BURP header information is saved as a string in the form LAT–LON–YYYYMMDDHHMMSTNID—-. Intention of this function is to be used for setting the unique identifier ‘code’ in struct_oss_obsdata. Can be called under the interface oss_obsdata_get_header_code.

Arguments

Arguments
  • lon [real ,in] :: longitude real (radians)

  • lat [real ,in] :: latitude real (radians)

  • date [integer ,in] :: date in YYYYMMDD

  • time [integer ,in] :: time in HHMM

  • stnid [character ,in] :: station ID

Return

code [character ] :: unique code

Call to

obsdata_get_header_code_i()

subroutine  obssubspacedata_mod/oss_obsdata_add_data1d(obsdata, val, code, maxsize[, dim1_opt])
Purpose

Add data value(s) to obsdata%data1d with associated identifier code.

Arguments

Arguments
  • obsdata [struct_oss_obsdata ,inout] :: Updated obsdata

  • val (*) [real ,in] :: data array to store in obsdata%data1d

  • code [character ,in] :: identifying code based on (lat,long,date,hhmm) if not also stnid

  • maxsize [integer ,in] :: max allowed size for obsdata dim1 value() dimension (optional)

  • obsdata :: Updated obsdata

Comments
  • Retrieval of values from obsdata%data1d to be done via oss_obsdata_get_element (or oss_obsdata_get_array1d).

  • If obsdata allocation is required for all processors (such as for use later with obsdata_MPIGather), allocation and/or initialization of arrays needs to be done at a corresponding appropriate location outside the obs operations such as in oss_setup to ensure allocation is done for all processors, including those without associated data. This is to ensure that rpn_comm_allgather will work in routine obsdata_MPIGather.

Options

dim1_opt [integer ,in,]

Call to

oss_obsdata_alloc(), utl_abort()

subroutine  obssubspacedata_mod/oss_obsdata_mpiallgather(obsdata)
Purpose

Gathers previously saved obsdata from all processors.

Arguments

Arguments

obsdata [struct_oss_obsdata ,inout] :: Local struct_oss_obsdata to become global

Comments
  • Assumes obsdata%dim1 (and obsdata%dim2) the same over all processors.

Called from

obsf_obssub_read()

Call to

mmpi_allgather_string(), oss_obsdata_dealloc(), oss_obsdata_alloc()

subroutine  obssubspacedata_mod/oss_get_comboidlist(obsspacedata, stnid_list, varno_list, unilev_list, num_elements, nset)
Purpose

Uses the subroutine oss_comboIdlist to compile a unique list of stnid, (stnid,varno) or (stnid,varno,multi/uni-level) combinations to be used in searches.

Arguments
obsSpaceData

Observation space data

Arguments
  • stnid_list (100) [character ,out] :: List of unique stnids

  • varno_list (100) [integer ,out] :: List of unique varno

  • unilev_list (100) [logical ,out] :: List of unique uni/multi-level identifications

  • num_elements [integer ,out] :: Number of unique elements in *_list arrrays

  • nset [integer ,out] ::

    Integer indicating grouping, with values indicating

    • 1: group by stnid

    • 2: group by (stnid,bufr)

    • 3: group by (stnid,bufr,multi/uni-level)

  • obsspacedata [struct_obs ,inout]

Called from

osd_obsdiagnostics()

Call to

oss_comboidlist(), obs_getheaderindex(), obs_headelem_i(), obs_bodyelem_i(), obs_elem_c()

subroutine  obssubspacedata_mod/oss_comboidlist([stnid_add_opt[, varno_add_opt[, unilev_add_opt[, stnid_list_opt[, varno_list_opt[, unilev_list_opt[, num_elements_opt[, initialize_opt[, nset_opt[, gather_mpi_opt[, all_combos_opt]]]]]]]]]]])
Purpose

Provide list of fixed or accumulated stnid, (stnid,varno) or (stnid,varno,multi/uni-level) combinations to be used in searches.

Can be used for both single processor and MPI mode, where ‘gather_mpi’ must be set to .true. at some point for use with MPI.

Called from osd_chem_diagnmostics in file obspacediag_mod.ftn90.

Arguments

Options
  • stnid_add_opt [character ,in,] :: stnid to add to stnid_list if part of unique set

  • varno_add_opt [integer ,in,] :: varno to add to varno_list if part of unique set

  • unilev_add_opt [logical ,in,] :: unilev logical to add to unilev_list if part of unique set

  • initialize_opt [logical ,in,] :: Initialize internal arrays and counters

  • gather_mpi_opt [logical ,in,] :: Gather all local MPI process and recompile unique lists

  • nset_opt [integer ,inout,] ::

    Integer indicating grouping of diagnostics. Input variable

    if initialize=.true., output variable otherwise. Values indicate

    • 1: group by stnid

    • 2: group by (stnid,bufr)

    • 3: group by (stnid,bufr,multi/uni-level)

    all combos_opt

    Indicates if all combinations specified by nset are to be use, or only those specified in the namelist NAMCHEM Input variable if initialize=.true., output variable otherwise.

  • stnid_list_opt (100) [character ,out,] :: List of unique stnids

  • varno_list_opt (100) [integer ,out,] :: List of unique varno

  • unilev_list_opt (100) [logical ,out,] :: List of unique uni/multi-level identifications

  • num_elements_opt [integer ,out,] :: Number of unique elements in *_list arrrays

  • all_combos_opt [logical ,inout,]

Called from

osd_obspostproc(), oss_get_comboidlist()

Call to

utl_abort(), utl_stnid_equal(), mmpi_allgather_string()

function  obssubspacedata_mod/oss_obsdata_code_len()
Purpose

Pass on oss_code_len value.

Return

oss_obsdata_code_len [integer ]

Called from

brpf_obssub_update()