Multi_Tool: Feature Relative by Lead using Multiple User-Defined Fields

model_applications/medium_range/ TCStat_SeriesAnalysis_fcstGFS _obsGFS_FeatureRelative _SeriesByLead_PyEmbed_Multiple_Diagnostics.conf

Scientific Objective

This use case calls multiple tools to produce diagnostic plots of systematic erros relative to a feature (e.g. hurricane, MCS, etc…). This use case calls two user provided python scripts that calculate diagnostics of interest (e.g. integrated vapor transport, potential vorticity, etc…). These user diagnostics are then used to define the systematic errors. This example calculates statistics over varying forecast leads with the ability to define lead groupings. This use case is very similar to the Multi_Tools: Feature Relative by Lead use case and the Multi_Tools: Feature Relative by Lead using User-Defined Fields. (ADeck,GFS:BDeck,GFS:ATCF,Grib2)

By maintaining focus of each evaluation time (or evaluation time series, in this case) on a user-defined area around a cyclone, the model statistical errors associated with cyclonic physical features (moisture flux, stability, strength of upper-level PV anomaly and jet, etc.) can be related directly to the model forecasts and provide improvement guidance by accurately depicting interactions with significant weather features around and within the cyclone. This is in contrast to the traditional method of regional averaging cyclone observations in a fixed grid, which “smooths out” system features and limits the meaningful metrics that can be gathered. Specifically, this use case creates bins of forecast lead times as specified by the given ranges which provides additional insight directly into forecast lead time accuracy.

Additionally, the ability to calculate model statistical errors based on user provided diagnostics allows the user to customize the feature relative analysis to suit their needs.

Datasets

This use case compares the Global Forecast System (GFS) forecast to the GFS analysis for hurricane Dorian. It is based on three user provided python scripts that calculate the diagnostic integrated vaport transport (IVT) baroclinic potential vorticity (PV), and saturation equivalent potential temperature (SEPT), respectively.

  • Variables required to calculate IVT: Levels required: all pressure levels >= 100mb #. Temperature #. v- component of wind #. u- component of wind #. Geopotential height #. Specific humidity OR Relative Humidity

  • Variables required to calculate PV: Levels required: all pressure levels >= 100mb #. U-wind #. V-wind #. Temperature

  • Variables required to calculate saturation equivalent potential temperature: Levels required: all pressure levels >= 100mb #. Temperature

  • Forecast dataset: GFS Grid 4 Forecast GFS Forecast data can be found at the following website: https://www.ncdc.noaa.gov/data-access/model-data/model-datasets/global-forcast-system-gfs - Initialization date: 20190830 - Initialization hours: 00, 06, 12, 18 UTC - Lead times: 90, 96, 102, 108, 114 - Format: Grib2 - Resolution: 0.5 degree

  • Observation dataset: GFS Grid 4 Analysis GFS Analysis data can be found at the following website: https://www.ncdc.noaa.gov/data-access/model-data/model-datasets/global-forcast-system-gfs - Valid date/time range: 20190902_18 - 20190904_12 every 6 hours - Format: Grib2 - Resolution: 0.5 degree

  • Hurricane Track Data Hurricane track data can be found at the following website: http://hurricanes.ral.ucar.edu/repository/data/ - ADeck Track File: aal052019.dat - BDeck Track File: bal052019.dat

External Dependencies

You will need to use a version of Python 3.7+ that has the following packages installed:

  • netCDF4

  • pygrib

  • cfgrib

  • metpy

  • xarray

If the version of Python used to compile MET did not have these libraries at the time of compilation, you will need to add these packages or create a new Python environment with these packages.

If this is the case, you will need to set the MET_PYTHON_EXE environment variable to the path of the version of Python you want to use. If you want this version of Python to only apply to this use case, set it in the [user_env_vars] section of a METplus configuration file.:

[user_env_vars]
MET_PYTHON_EXE = /path/to/python/with/required/packages/bin/python

METplus Components

This use case first runs PyEmbedIngest to run the user provided python scripts to calculate the desired diagnostics (in this example, IVT, PV and SEPT). PyEmbedIngest runs the RegridDataPlane tool to write IVT, PV, and SEPTto a MET readable netCDF file. Then TCPairs and ExtractTiles are run to generate matched tropical cyclone data and regrid them into appropriately-sized tiles along a storm track. The MET tc-stat tool is used to filter the track data and the MET regrid-dataplane tool is used to regrid the data (GRIB1 or GRIB2 into netCDF). Next, a series analysis by lead time is performed on the results and plots (.ps and .png) are generated for all variable-level-stat combinations from the specified variables, levels, and requested statistics. If lead grouping is turned on, the final results are aggregated into forecast hour groupings as specified by the start, end and increment in the METplus configuration file, as well as labels to identify each forecast hour grouping. If lead grouping is not turned out the final results will be written out for each requested lead time.

METplus Workflow

This use case loops by process which means that each tool is run for all times before moving to the next tool. The tool order is as follows:

PyEmbedIngest, TCPairs, ExtractTiles, SeriesByLead

This example loops by forecast/lead time (with begin, end, and increment as specified in the METplus TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_Multiple_Diagnostics.conf file).

4 initialization times will be run over 5 lead times:

Init: 20190830_00Z
Forecast lead: 90, 96, 102, 108, 114

Init: 20190830_06Z
Forecast lead: 90, 96, 102, 108, 114

Init: 20190830_12Z
Forecast lead: 90, 96, 102, 108, 114

Init: 20190830_18Z
Forecast lead: 90, 96, 102, 108, 114

METplus Configuration

METplus first loads all of the configuration files found in parm/metplus_config, then it loads any configuration files passed to METplus via the command line with the -c option, i.e. -c parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_Multiple_Diagnostics.conf

[config]

# Documentation for this use case can be found at
# https://metplus.readthedocs.io/en/latest/generated/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics.html

# For additional information, please see the METplus Users Guide.
# https://metplus.readthedocs.io/en/latest/Users_Guide

###
# Processes to run
# https://metplus.readthedocs.io/en/latest/Users_Guide/systemconfiguration.html#process-list
###

PROCESS_LIST = PyEmbedIngest, TCPairs, TCStat, ExtractTiles, TCStat(for_series_analysis), SeriesAnalysis


###
# Time Info
# LOOP_BY options are INIT, VALID, RETRO, and REALTIME
# If set to INIT or RETRO:
#   INIT_TIME_FMT, INIT_BEG, INIT_END, and INIT_INCREMENT must also be set
# If set to VALID or REALTIME:
#   VALID_TIME_FMT, VALID_BEG, VALID_END, and VALID_INCREMENT must also be set
# LEAD_SEQ is the list of forecast leads to process
# https://metplus.readthedocs.io/en/latest/Users_Guide/systemconfiguration.html#timing-control
###

LOOP_BY = INIT
INIT_TIME_FMT = %Y%m%d%H
INIT_BEG = 2019083000
INIT_END = 2019083023
INIT_INCREMENT = 21600

LEAD_SEQ = 90, 96, 102, 108, 114

TC_PAIRS_RUNTIME_FREQ = RUN_ONCE
SERIES_ANALYSIS_RUNTIME_FREQ = RUN_ONCE_PER_LEAD

SERIES_ANALYSIS_RUN_ONCE_PER_STORM_ID = False


###
# File I/O
# https://metplus.readthedocs.io/en/latest/Users_Guide/systemconfiguration.html#directory-and-filename-template-info
###

PY_EMBED_INGEST_1_OUTPUT_DIR = {OUTPUT_BASE}/py_embed_out
PY_EMBED_INGEST_1_OUTPUT_TEMPLATE = {init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.nc

PY_EMBED_INGEST_2_OUTPUT_DIR = {PY_EMBED_INGEST_1_OUTPUT_DIR}
PY_EMBED_INGEST_2_OUTPUT_TEMPLATE = {valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.nc


TC_PAIRS_ADECK_INPUT_DIR = {INPUT_BASE}/model_applications/medium_range/dorian_data/track_data
TC_PAIRS_ADECK_TEMPLATE = a{basin?fmt=%s}052019.dat

TC_PAIRS_BDECK_INPUT_DIR = {TC_PAIRS_ADECK_INPUT_DIR}
TC_PAIRS_BDECK_TEMPLATE = b{basin?fmt=%s}052019.dat

TC_PAIRS_REFORMAT_DIR = {OUTPUT_BASE}/track_data_atcf
TC_PAIRS_SKIP_IF_REFORMAT_EXISTS = no

TC_PAIRS_OUTPUT_DIR = {OUTPUT_BASE}/tc_pairs
TC_PAIRS_OUTPUT_TEMPLATE = {basin?fmt=%s}q{INIT_BEG}.dorian
TC_PAIRS_SKIP_IF_OUTPUT_EXISTS = no


TC_STAT_LOOKIN_DIR = {TC_PAIRS_OUTPUT_DIR}

TC_STAT_OUTPUT_DIR = {EXTRACT_TILES_OUTPUT_DIR}
TC_STAT_DUMP_ROW_TEMPLATE = filter_{init?fmt=%Y%m%d_%H}.tcst


EXTRACT_TILES_TC_STAT_INPUT_DIR = {TC_STAT_OUTPUT_DIR}
EXTRACT_TILES_TC_STAT_INPUT_TEMPLATE = {TC_STAT_DUMP_ROW_TEMPLATE}

EXTRACT_TILES_GRID_INPUT_DIR = {PY_EMBED_INGEST_1_OUTPUT_DIR}
FCST_EXTRACT_TILES_INPUT_DIR = {PY_EMBED_INGEST_1_OUTPUT_DIR}
FCST_EXTRACT_TILES_INPUT_TEMPLATE = {init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.nc

OBS_EXTRACT_TILES_INPUT_DIR = {PY_EMBED_INGEST_1_OUTPUT_DIR}
OBS_EXTRACT_TILES_INPUT_TEMPLATE = {valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.nc

EXTRACT_TILES_OUTPUT_DIR = {OUTPUT_BASE}/extract_tiles
FCST_EXTRACT_TILES_OUTPUT_TEMPLATE = {init?fmt=%Y%m%d_%H}/{storm_id}/FCST_TILE_F{lead?fmt=%3H}_{MODEL}_gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.nc
OBS_EXTRACT_TILES_OUTPUT_TEMPLATE = {init?fmt=%Y%m%d_%H}/{storm_id}/OBS_TILE_F{lead?fmt=%3H}_{MODEL}_gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.nc
EXTRACT_TILES_SKIP_IF_OUTPUT_EXISTS = no


FCST_SERIES_ANALYSIS_INPUT_DIR = {EXTRACT_TILES_OUTPUT_DIR}
FCST_SERIES_ANALYSIS_INPUT_TEMPLATE = {FCST_EXTRACT_TILES_OUTPUT_TEMPLATE}

OBS_SERIES_ANALYSIS_INPUT_DIR = {EXTRACT_TILES_OUTPUT_DIR}
OBS_SERIES_ANALYSIS_INPUT_TEMPLATE = {OBS_EXTRACT_TILES_OUTPUT_TEMPLATE}

SERIES_ANALYSIS_TC_STAT_INPUT_DIR = {SERIES_ANALYSIS_OUTPUT_DIR}
SERIES_ANALYSIS_TC_STAT_INPUT_TEMPLATE = {TC_STAT_DUMP_ROW_TEMPLATE}

SERIES_ANALYSIS_OUTPUT_DIR = {OUTPUT_BASE}/series_analysis_lead
SERIES_ANALYSIS_OUTPUT_TEMPLATE = {label}/series_F{fcst_beg}_{fcst_name}_{fcst_level}.nc


###
# Field Info
# https://metplus.readthedocs.io/en/latest/Users_Guide/systemconfiguration.html#field-info
###

# Used by extract tiles and series analysis to define the records of
#  interest to be retrieved from the grib2 file 

MODEL = GFSO

BOTH_VAR1_NAME = ivt
BOTH_VAR1_LEVELS = Surface

BOTH_VAR2_NAME = pv
BOTH_VAR2_LEVELS = Surface

BOTH_VAR3_NAME = sept
BOTH_VAR3_LEVELS = Surface


###
# PyEmbedIngest Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#pyembedingest
###

CONFIG_DIR={PARM_BASE}/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics
MODEL_DIR = {INPUT_BASE}/model_applications/medium_range/dorian_data/model_data

# 1st INGEST INSTANCE: Forecast

# IVT
PY_EMBED_INGEST_1_SCRIPT_1 = {CONFIG_DIR}/gfs_ivt_fcst.py {MODEL_DIR}/{init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.grb2
PY_EMBED_INGEST_1_OUTPUT_FIELD_NAME_1 = ivt

# PV
PY_EMBED_INGEST_1_SCRIPT_2 = {CONFIG_DIR}/gfs_pv_fcst.py {MODEL_DIR}/{init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.grb2
PY_EMBED_INGEST_1_OUTPUT_FIELD_NAME_2 = pv

# SEPT
PY_EMBED_INGEST_1_SCRIPT_3 = {CONFIG_DIR}/gfs_sept_fcst.py {MODEL_DIR}/{init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.grb2
PY_EMBED_INGEST_1_OUTPUT_FIELD_NAME_3 = sept

PY_EMBED_INGEST_1_TYPE = NUMPY

PY_EMBED_INGEST_1_OUTPUT_GRID = {MODEL_DIR}/{init?fmt=%Y%m%d}/gfs_4_{init?fmt=%Y%m%d}_{init?fmt=%H}00_{lead?fmt=%3H}.grb2

# 2nd INGEST INSTANCE: Analysis

# IVT
PY_EMBED_INGEST_2_SCRIPT_1 = {CONFIG_DIR}/gfs_ivt_analysis.py {MODEL_DIR}/{valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.grb2
PY_EMBED_INGEST_2_OUTPUT_FIELD_NAME_1 = ivt

# PV
PY_EMBED_INGEST_2_SCRIPT_2 = {CONFIG_DIR}/gfs_pv_analysis.py {MODEL_DIR}/{valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.grb2
PY_EMBED_INGEST_2_OUTPUT_FIELD_NAME_2 = pv

# SEPT
PY_EMBED_INGEST_2_SCRIPT_3 = {CONFIG_DIR}/gfs_sept_analysis.py {MODEL_DIR}/{valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.grb2
PY_EMBED_INGEST_2_OUTPUT_FIELD_NAME_3 = sept

PY_EMBED_INGEST_2_TYPE = NUMPY

PY_EMBED_INGEST_2_OUTPUT_GRID = {MODEL_DIR}/{valid?fmt=%Y%m%d}/gfs_4_{valid?fmt=%Y%m%d}_{valid?fmt=%H}00_000.grb2


###
# TCPairs Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#tcpairs
###

TC_PAIRS_INIT_INCLUDE =
TC_PAIRS_INIT_EXCLUDE =

TC_PAIRS_VALID_BEG =
TC_PAIRS_VALID_END =

TC_PAIRS_STORM_ID =
TC_PAIRS_BASIN =
TC_PAIRS_CYCLONE =
TC_PAIRS_STORM_NAME =

TC_PAIRS_DLAND_FILE = {MET_INSTALL_DIR}/share/met/tc_data/dland_global_tenth_degree.nc

TC_PAIRS_REFORMAT_DECK = no
TC_PAIRS_REFORMAT_TYPE = SBU

TC_PAIRS_MISSING_VAL_TO_REPLACE = -99
TC_PAIRS_MISSING_VAL = -9999


###
# TCStat Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#tcstat
###

TC_STAT_JOB_ARGS = -job filter -basin AL -dump_row {TC_STAT_OUTPUT_DIR}/{TC_STAT_DUMP_ROW_TEMPLATE}

TC_STAT_MATCH_POINTS = true

TC_STAT_AMODEL = {MODEL}
TC_STAT_BMODEL =
TC_STAT_DESC =
TC_STAT_STORM_ID =
TC_STAT_BASIN =
TC_STAT_CYCLONE =
TC_STAT_STORM_NAME =

TC_STAT_INIT_BEG =
TC_STAT_INIT_END =
TC_STAT_INIT_INCLUDE = {init?fmt=%Y%m%d_%H}
TC_STAT_INIT_EXCLUDE =
TC_STAT_INIT_HOUR =

TC_STAT_VALID_BEG =
TC_STAT_VALID_END = 
TC_STAT_VALID_INCLUDE =
TC_STAT_VALID_EXCLUDE =
TC_STAT_VALID_HOUR =
TC_STAT_LEAD_REQ =
TC_STAT_INIT_MASK =
TC_STAT_VALID_MASK =

TC_STAT_VALID_HOUR =
TC_STAT_LEAD =

TC_STAT_TRACK_WATCH_WARN =

TC_STAT_COLUMN_THRESH_NAME =
TC_STAT_COLUMN_THRESH_VAL =

TC_STAT_COLUMN_STR_NAME =
TC_STAT_COLUMN_STR_VAL =

TC_STAT_INIT_THRESH_NAME =
TC_STAT_INIT_THRESH_VAL =

TC_STAT_INIT_STR_NAME =
TC_STAT_INIT_STR_VAL =

TC_STAT_WATER_ONLY =

TC_STAT_LANDFALL =

TC_STAT_LANDFALL_BEG =
TC_STAT_LANDFALL_END =


###
# ExtractTiles Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#extracttiles
###

EXTRACT_TILES_NLAT = 60
EXTRACT_TILES_NLON = 60

EXTRACT_TILES_DLAT = 0.5
EXTRACT_TILES_DLON = 0.5

EXTRACT_TILES_LON_ADJ = 15
EXTRACT_TILES_LAT_ADJ = 15


###
# TCStat (for SeriesAnalysis) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#tcstat
###

# Settings specific to the TCStat(for_series_analysis) process that was set
# in the PROCESS_LIST. Any TC_STAT_* variable not set in this section will use
# the value set outside of this section

[for_series_analysis]

TC_STAT_JOB_ARGS = -job filter -init_beg {INIT_BEG} -init_end {INIT_END} -dump_row {TC_STAT_OUTPUT_DIR}/{TC_STAT_DUMP_ROW_TEMPLATE}

TC_STAT_OUTPUT_DIR = {SERIES_ANALYSIS_OUTPUT_DIR}
TC_STAT_LOOKIN_DIR = {EXTRACT_TILES_OUTPUT_DIR}


###
# SeriesAnalysis Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#seriesanalysis
###

[config]

SERIES_ANALYSIS_REGRID_TO_GRID = FCST
SERIES_ANALYSIS_REGRID_METHOD = FORCE

SERIES_ANALYSIS_STAT_LIST = TOTAL, FBAR, OBAR, ME

SERIES_ANALYSIS_BACKGROUND_MAP = yes

SERIES_ANALYSIS_BLOCK_SIZE = 4000

SERIES_ANALYSIS_IS_PAIRED = True

SERIES_ANALYSIS_GENERATE_PLOTS = yes
SERIES_ANALYSIS_GENERATE_ANIMATIONS = yes

PLOT_DATA_PLANE_TITLE = {MODEL} series_F{fcst_beg} Forecasts{nseries}, {stat} for {fcst_name} {fcst_level}


[user_env_vars]

PV_LAYER_MIN_PRESSURE=100.0
PV_LAYER_MAX_PRESSURE=1000.0
IVT_LAYER_MIN_PRESSURE=100.0
IVT_LAYER_MAX_PRESSURE=1000.0
SEPT_LAYER_MIN_PRESSURE=100.0
SEPT_LAYER_MAX_PRESSURE=1000.0

MET Configuration

METplus sets environment variables based on user settings in the METplus configuration file. See How METplus controls MET config file settings for more details.

YOU SHOULD NOT SET ANY OF THESE ENVIRONMENT VARIABLES YOURSELF! THEY WILL BE OVERWRITTEN BY METPLUS WHEN IT CALLS THE MET TOOLS!

If there is a setting in the MET configuration file that is currently not supported by METplus you’d like to control, please refer to: Overriding Unsupported MET config file settings

TCPairsConfig_wrapped

Note

See the TCPairs MET Configuration section of the User’s Guide for more information on the environment variables used in the file below:

////////////////////////////////////////////////////////////////////////////////
//
// Default TCPairs configuration file
//
////////////////////////////////////////////////////////////////////////////////

//
// ATCF file format reference:
//   http://www.nrlmry.navy.mil/atcf_web/docs/database/new/abrdeck.html
//

//
// Models
//
//model = [
${METPLUS_MODEL}

//
// Description
//
//desc =
${METPLUS_DESC}

//
// Storm identifiers
//
//storm_id = [
${METPLUS_STORM_ID}

//
// Basins
//
//basin = [
${METPLUS_BASIN}

//
// Cyclone numbers
//
//cyclone = [
${METPLUS_CYCLONE}

//
// Storm names
//
//storm_name = [
${METPLUS_STORM_NAME}

//
// Model initialization time windows to include or exclude
//
//init_beg =
${METPLUS_INIT_BEG}
//init_end =
${METPLUS_INIT_END}
// init_inc = [
${METPLUS_INIT_INC}
// init_exc = [
${METPLUS_INIT_EXC}

//
// Valid model time windows to include or exclude
//
//valid_beg =
${METPLUS_VALID_BEG}
//valid_end =
${METPLUS_VALID_END}
// valid_inc = [
${METPLUS_VALID_INC}
// valid_exc = [
${METPLUS_VALID_EXC}

// write_valid =
${METPLUS_WRITE_VALID}

//
// Model initialization hours
//
init_hour = [];

//
// Required lead time in hours
//
lead_req = [];

//
// lat/lon polylines defining masking regions
//
init_mask  = "";
valid_mask = "";

//
// Specify if the code should check for duplicate ATCF lines
//
//check_dup =
${METPLUS_CHECK_DUP}


//
// Specify special processing to be performed for interpolated models.
// Set to NONE, FILL, or REPLACE.
//
//interp12 =
${METPLUS_INTERP12}

//
// Specify how consensus forecasts should be defined
//
//consensus =
${METPLUS_CONSENSUS_LIST}


//
// Forecast lag times
//
lag_time = [];

//
// CLIPER/SHIFOR baseline forecasts to be derived from the BEST
// and operational (CARQ) tracks.
//
best_technique = [ "BEST" ];
best_baseline  = [];
oper_technique = [ "CARQ" ];
oper_baseline  = [];

//
// Specify the datasets to be searched for analysis tracks (NONE, ADECK, BDECK,
// or BOTH).
//
anly_track = BDECK;

//
// Specify if only those track points common to both the ADECK and BDECK
// tracks be written out.
//
//match_points =
${METPLUS_MATCH_POINTS}

//
// Specify the NetCDF output of the gen_dland tool containing a gridded
// representation of the minimum distance to land.
//
//dland_file =
${METPLUS_DLAND_FILE}

//
// Specify watch/warning information:
//   - Input watch/warning filename
//   - Watch/warning time offset in seconds
//
watch_warn = {
   file_name   = "MET_BASE/tc_data/wwpts_us.txt";
   time_offset = -14400;
}


//diag_info_map = {
${METPLUS_DIAG_INFO_MAP_LIST}

//diag_convert_map = {
${METPLUS_DIAG_CONVERT_MAP_LIST}

//
// Indicate a version number for the contents of this configuration file.
// The value should generally not be modified.
//
//version = "V9.0";

tmp_dir = "${MET_TMP_DIR}";

${METPLUS_MET_CONFIG_OVERRIDES}

TCStatConfig_wrapped

Note

See the TCStat MET Configuration section of the User’s Guide for more information on the environment variables used in the file below:

///////////////////////////////////////////////////////////////////////////////
//
// Default TCStat configuration file
//
////////////////////////////////////////////////////////////////////////////////

//
// The parameters listed below are used to filter the TC-STAT data down to the
// desired subset of lines over which statistics are to be computed.  Only
// those lines which meet ALL of the criteria specified will be retained.
//
// The settings that are common to all jobs may be specified once at the top
// level.  If no selection is listed for a parameter, that parameter will not
// be used for filtering.  If multiple selections are listed for a parameter,
// the analyses will be performed on their union.
//

//
// Stratify by the AMODEL or BMODEL columns.
//
//amodel = [
${METPLUS_AMODEL}
//bmodel = [
${METPLUS_BMODEL}

//
// Stratify by the DESC column.
//
//desc = [
${METPLUS_DESC}

//
// Stratify by the STORM_ID column.
//
//storm_id = [
${METPLUS_STORM_ID}

//
// Stratify by the BASIN column.
// May add using the "-basin" job command option.
//
//basin = [
${METPLUS_BASIN}

//
// Stratify by the CYCLONE column.
// May add using the "-cyclone" job command option.
//
//cyclone = [
${METPLUS_CYCLONE}

//
// Stratify by the STORM_NAME column.
// May add using the "-storm_name" job command option.
//
//storm_name = [
${METPLUS_STORM_NAME}

//
// Stratify by the INIT times.
// Model initialization time windows to include or exclude
// May modify using the "-init_beg", "-init_end", "-init_inc",
// and "-init_exc" job command options.
//
//init_beg =
${METPLUS_INIT_BEG}
//init_end =
${METPLUS_INIT_END}
//init_inc = [
${METPLUS_INIT_INC}
//init_exc = [
${METPLUS_INIT_EXC}

//
// Stratify by the VALID times.
//
//valid_beg =
${METPLUS_VALID_BEG}
//valid_end =
${METPLUS_VALID_END}
//valid_inc = [
${METPLUS_VALID_INC}
//valid_exc = [
${METPLUS_VALID_EXC}

//
// Stratify by the initialization and valid hours and lead time.
//
//init_hour = [
${METPLUS_INIT_HOUR}
//valid_hour = [
${METPLUS_VALID_HOUR}
//lead = [
${METPLUS_LEAD}

//
// Select tracks which contain all required lead times.
//
//lead_req = [
${METPLUS_LEAD_REQ}

//
// Stratify by the INIT_MASK and VALID_MASK columns.
//
//init_mask = [
${METPLUS_INIT_MASK}
//valid_mask = [
${METPLUS_VALID_MASK}


//
// Stratify by the LINE_TYPE column.
//
//line_type =
${METPLUS_LINE_TYPE}


//
// Stratify by checking the watch/warning status for each track point
// common to both the ADECK and BDECK tracks.  If the watch/warning status
// of any of the track points appears in the list, retain the entire track.
//
//track_watch_warn = [
${METPLUS_TRACK_WATCH_WARN}

//
// Stratify by applying thresholds to numeric data columns.
//
//column_thresh_name = [
${METPLUS_COLUMN_THRESH_NAME}
//column_thresh_val = [
${METPLUS_COLUMN_THRESH_VAL}

//
// Stratify by performing string matching on non-numeric data columns.
//
//column_str_name = [
${METPLUS_COLUMN_STR_NAME}
//column_str_val = [
${METPLUS_COLUMN_STR_VAL}

//
// Stratify by excluding strings in non-numeric data columns.
//
//column_str_exc_name =
${METPLUS_COLUMN_STR_EXC_NAME}

//column_str_exc_val =
${METPLUS_COLUMN_STR_EXC_VAL}

//
// Similar to the column_thresh options above
//
//init_thresh_name = [
${METPLUS_INIT_THRESH_NAME}
//init_thresh_val = [
${METPLUS_INIT_THRESH_VAL}

//
// Similar to the column_str options above
//
//init_str_name = [
${METPLUS_INIT_STR_NAME}
//init_str_val = [
${METPLUS_INIT_STR_VAL}

//
// Similar to the column_str_exc options above
//
//init_str_exc_name =
${METPLUS_INIT_STR_EXC_NAME}

//init_str_exc_val =
${METPLUS_INIT_STR_EXC_VAL}

//diag_thresh_name =
${METPLUS_DIAG_THRESH_NAME}

//diag_thresh_val =
${METPLUS_DIAG_THRESH_VAL}

//init_diag_thresh_name =
${METPLUS_INIT_DIAG_THRESH_NAME}

//init_diag_thresh_val =
${METPLUS_INIT_DIAG_THRESH_VAL}


//
// Stratify by the ADECK and BDECK distances to land.
//
//water_only =
${METPLUS_WATER_ONLY}

//
// Specify whether only those track points occurring near landfall should be
// retained, and define the landfall retention window in HH[MMSS] format
// around the landfall time.
//
//landfall =
${METPLUS_LANDFALL}
//landfall_beg =
${METPLUS_LANDFALL_BEG}
//landfall_end =
${METPLUS_LANDFALL_END}

//
// Specify whether only those track points common to both the ADECK and BDECK
// tracks should be retained.  May modify using the "-match_points" job command
// option.
//
//match_points =
${METPLUS_MATCH_POINTS}

//event_equal =
${METPLUS_EVENT_EQUAL}

//event_equal_lead =
${METPLUS_EVENT_EQUAL_LEAD}

//out_init_mask =
${METPLUS_OUT_INIT_MASK}

//out_valid_mask =
${METPLUS_OUT_VALID_MASK}


//
// Array of TCStat analysis jobs to be performed on the filtered data
//
//jobs = [
${METPLUS_JOBS}

tmp_dir = "${MET_TMP_DIR}";

${METPLUS_MET_CONFIG_OVERRIDES}

SeriesAnalysisConfig_wrapped

Note

See the SeriesAnalysis MET Configuration section of the User’s Guide for more information on the environment variables used in the file below:

////////////////////////////////////////////////////////////////////////////////
//
// Series-Analysis configuration file.
//
// For additional information, see the MET_BASE/config/README file.
//
////////////////////////////////////////////////////////////////////////////////

//
// Output model name to be written
//
//model =
${METPLUS_MODEL}

//
// Output description to be written
//
//desc =
${METPLUS_DESC}

//
// Output observation type to be written
//
//obtype =
${METPLUS_OBTYPE}

////////////////////////////////////////////////////////////////////////////////

//
// Verification grid
// May be set separately in each "field" entry
//
//regrid = {
${METPLUS_REGRID_DICT}

////////////////////////////////////////////////////////////////////////////////

censor_thresh = [];
censor_val    = [];
//cat_thresh =
${METPLUS_CAT_THRESH}
cnt_thresh    = [ NA ];
cnt_logic     = UNION;

//
// Forecast and observation fields to be verified
//
fcst = {
   ${METPLUS_FCST_FILE_TYPE}
   ${METPLUS_FCST_CAT_THRESH}
   //field = [
   ${METPLUS_FCST_FIELD}
}
obs = {
   ${METPLUS_OBS_FILE_TYPE}
   ${METPLUS_OBS_CAT_THRESH}
   //field = [
   ${METPLUS_OBS_FIELD}
}

////////////////////////////////////////////////////////////////////////////////

//
// Climatology data
//
//climo_mean = {
${METPLUS_CLIMO_MEAN_DICT}


//climo_stdev = {
${METPLUS_CLIMO_STDEV_DICT}

//climo_cdf = {
${METPLUS_CLIMO_CDF_DICT}

////////////////////////////////////////////////////////////////////////////////

//
// Confidence interval settings
//
ci_alpha  = [ 0.05 ];

boot = {
   interval = PCTILE;
   rep_prop = 1.0;
   n_rep    = 0;
   rng      = "mt19937";
   seed     = "";
}

////////////////////////////////////////////////////////////////////////////////

//
// Verification masking regions
//
//mask = {
${METPLUS_MASK_DICT}

//
// Number of grid points to be processed concurrently.  Set smaller to use
// less memory but increase the number of passes through the data.
//
//block_size =
${METPLUS_BLOCK_SIZE}

//
// Ratio of valid matched pairs to compute statistics for a grid point
//
//vld_thresh =
${METPLUS_VLD_THRESH}

////////////////////////////////////////////////////////////////////////////////

//
// Statistical output types
//
//output_stats = {
${METPLUS_OUTPUT_STATS_DICT}

////////////////////////////////////////////////////////////////////////////////

//hss_ec_value =
${METPLUS_HSS_EC_VALUE}
rank_corr_flag = FALSE;

tmp_dir = "${MET_TMP_DIR}";

//version        = "V10.0";

////////////////////////////////////////////////////////////////////////////////

${METPLUS_MET_CONFIG_OVERRIDES}

Python Embedding

This use case uses four Python embedding scripts to read input data, two for the forecast data and two for the analysis data. The multiple datatype input requires the two-script approach.

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_ivt_fcst.py

# This script is a combination of two scripts originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR. 
# May 2020
###################################################################################################

import pygrib
import numpy as np
import sys
import os
import re
import datetime as dt
import metpy.calc as mc

###################################################################################################

def ivt(input_file):
    grbs = pygrib.open(input_file)
    g = 9.81    # Setting gravity constant
    print(input_file)
    grbs.rewind()
    
    # Initialize variable arrays 
    levs = [] # Levels 
    q    = [] # Specific humidity
    hgt  = [] # Geopotential height
    temp = [] # Temperature
    u    = [] # u-wind
    v    = [] # v-wind

    # First obtain the levels we will use
    # These are in hPa in the file, so directly compare with user supplied min/max
    levs = sorted(set([grb.level for grb in grbs if float(grb.level) >= float(os.environ.get('IVT_LAYER_MIN_PRESSURE',100.0)) and float(grb.level) <= float(os.environ.get('IVT_LAYER_MAX_PRESSURE',1000.0))]))
    
    # Fill in variable arrays from input file.
    grbs.rewind()
    for grb in grbs:
        if not grb.level in levs:
            continue
        elif np.logical_and('v-' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            v.append(grb.values)        
        elif np.logical_and('u-' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            u.append(grb.values)
        elif np.logical_and('Temperature' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            temp.append(grb.values)
        elif np.logical_and('Geopotential' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            hgt.append(grb.values)
        elif np.logical_and('Specific' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            q.append(grb.values)
        
    temp = np.array(temp)
    hgt  = np.array(hgt)
    u    = np.array(u)
    v    = np.array(v)
    
    grbs.rewind()
    
    # If we didn't find specific humidity, look for relative humidity.
    if len(q) == 0:
        for grb in grbs:
            if not grb.level in levs:
                continue
            if np.logical_and('Relative' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
                q.append(grb.values)

        levs = np.array(levs)
        # Clausius-Clapeyron time
        es = 610.78*np.exp((17.67*(temp-273.15)/(temp-29.65))) # Calculate saturation vapor pressure
        e = es*(np.array(q)/100) # Calculate vapor pressure 
        w = 0.622*es/(levs[:,None,None]*100) # Calculate water vapor
        q = w/(w+1) # Calculate specific humidity
    q = np.array(q)
    
    uv = np.sqrt(u**2+v**2) # Calculate wind
    mflux_total = np.sum(q,axis=0)*(1/g)*np.mean(uv,axis=0)*(np.max(levs)-np.min(levs)) #calculate mass flux
    met_data = mflux_total.copy() #Pass mass flux to be used by MET tools
    print(np.max(met_data))
    #np.save('{}.npy'.format(sys.argv[1]),mflux_total)
    grbs.close()

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = ivt(input_file) #Call function to calculate IVT

met_data = data
met_data = met_data.astype('float64')

# Automatically fill out time information from input file.
file_regex = r"^.*([0-9]{8}_[0-9]{4})_([0-9]{3}).*$"
match = re.match(file_regex,
                 os.path.basename(input_file).replace('-', '_'))
if not match:
    print(f"Could not extract time information from filename: {input_file} using regex {file_regex}")
    sys.exit(1)

init = dt.datetime.strptime(match.group(1), '%Y%m%d_%H%M')
lead = int(match.group(2))
valid = init + dt.timedelta(hours=lead)

print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

print(valid)
print(init)
print(lead)


attrs = {
   'valid': valid.strftime("%Y%m%d_%H%M%S"),
   'init':  init.strftime("%Y%m%d_%H%M%S"),
   'lead':  str(int(lead)),
   'accum': '00',

   'name':      'ivt',
   'long_name': 'integrated_vapor_transport',
   'level':     'Surface',
   'units':     'UNKNOWN',

   'grid': {
       'name': 'Global 0.5 Degree',
       'type' :   'LatLon',
       'lat_ll' : -90.0,
       'lon_ll' : 0.0,
       'delta_lat' :   0.5,
       'delta_lon' :   0.5,
       'Nlat' :      361,
       'Nlon' :      720,
   }
}

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_pv_fcst.py

# This script calculates potential vorticity (PV) from variables found in the GFS forecast model grib files. This script is originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR.
# July 2020
###################################################################################################

import sys
import os
import re
import datetime as dt
from metpy import calc as mpcalc
from metpy.units import units
import xarray as xr
import cfgrib

###################################################################################################

def pv(input_file):

    # Vars
    grib_vars = ['t','u','v']

    # Load a list of datasets, one for each variable we want
    ds_list = [cfgrib.open_datasets(input_file,backend_kwargs={'filter_by_keys':{'typeOfLevel':'isobaricInhPa','shortName':v},'indexpath':''}) for v in grib_vars]

    # Flatten the list of lists to a single list of datasets
    ds_flat = [x.sel(isobaricInhPa=x.isobaricInhPa[x.isobaricInhPa>=100.0].values) for ds in ds_list for x in ds]

    # Merge the variables into a single dataset
    ds = xr.merge(ds_flat)

    # Add pressure
    ds['p'] = xr.DataArray(ds.isobaricInhPa.values,dims=['isobaricInhPa'],coords={'isobaricInhPa':ds.isobaricInhPa.values},attrs={'units':'hPa'}).broadcast_like(ds['t']) 

    # Calculate potential temperature
    ds['theta'] = mpcalc.potential_temperature(ds['p'].metpy.convert_units('Pa'),ds['t'])
  
    # Compute baroclinic PV
    ds['pv'] = mpcalc.potential_vorticity_baroclinic(ds['theta'],ds['p'].metpy.convert_units('Pa'),ds['u'],ds['v'],latitude=ds.latitude)/(1.0e-6)
    
    met_data = ds['pv'].sel(isobaricInhPa=slice(float(os.environ.get('PV_LAYER_MAX_PRESSURE',1000.0)),float(os.environ.get('PV_LAYER_MIN_PRESSURE',100.0)))).mean(axis=0).values

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = pv(input_file) #Call function to calculate PV

met_data = data
met_data = met_data.astype('float64')

print("max", data.max())
print("min", data.min())

# Automatically fill out time information from input file.
file_regex = r"^.*([0-9]{8}_[0-9]{4})_([0-9]{3}).*$"
match = re.match(file_regex, os.path.basename(input_file).replace('-', '_'))
if not match:
    print(f"Could not extract time information from filename: {input_file} using regex {file_regex}")
    sys.exit(1)
    
init = dt.datetime.strptime(match.group(1), '%Y%m%d_%H%M')
lead = int(match.group(2))
valid = init + dt.timedelta(hours=lead)

print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

print(valid)
print(init)
print(lead)

attrs = {
        'valid': valid.strftime("%Y%m%d_%H%M%S"),
        'init':  init.strftime("%Y%m%d_%H%M%S"),
        'lead':  str(int(lead)),
        'accum': '00',
        
        'name':      'pv',
        'long_name': 'potential_vorticity',
        'level':     'Surface',
        'units':     'PV Units',
        
        'grid': {
            'name': 'Global 0.5 Degree',
            'type' :   'LatLon',
            'lat_ll' : -90.0,
            'lon_ll' : 0.0,
            'delta_lat' :   0.5,
            'delta_lon' :   0.5,
            'Nlat' :      361,
            'Nlon' :      720,
            }
        }

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_sept_fcst.py

# This script calculates potential vorticity (PV) from variables found in the GFS forecast model grib files. This script is originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR.
# July 2020
###################################################################################################

import sys
import os
import re
import datetime as dt
from metpy import calc as mpcalc
from metpy.units import units
import xarray as xr
import cfgrib

###################################################################################################

def sept(input_file):

    # Vars
    grib_vars = ['t']

    # Load a list of datasets, one for each variable we want
    ds_list = [cfgrib.open_datasets(input_file,backend_kwargs={'filter_by_keys':{'typeOfLevel':'isobaricInhPa','shortName':v},'indexpath':''}) for v in grib_vars]

    # Flatten the list of lists to a single list of datasets
    ds_flat = [x.sel(isobaricInhPa=x.isobaricInhPa[x.isobaricInhPa>=100.0].values) for ds in ds_list for x in ds]

    # Merge the variables into a single dataset
    ds = xr.merge(ds_flat)

    # Add pressure
    ds['p'] = xr.DataArray(ds.isobaricInhPa.values,dims=['isobaricInhPa'],coords={'isobaricInhPa':ds.isobaricInhPa.values},attrs={'units':'hPa'}).broadcast_like(ds['t'])

    # Calculate saturation equivalent potential temperature
    ds['sept'] = mpcalc.saturation_equivalent_potential_temperature(ds['p'].metpy.convert_units('Pa'),ds['t'])

    met_data = ds['sept'].sel(isobaricInhPa=slice(float(os.environ.get('SEPT_LAYER_MAX_PRESSURE',1000.0)),float(os.environ.get('SEPT_LAYER_MIN_PRESSURE',100.0)))).mean(axis=0).values 

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = sept(input_file) #Call function to calculate PV

met_data = data
met_data = met_data.astype('float64')

print("max", data.max())
print("min", data.min())

# Automatically fill out time information from input file.
file_regex = r"^.*([0-9]{8}_[0-9]{4})_([0-9]{3}).*$"
match = re.match(file_regex, os.path.basename(input_file).replace('-', '_'))
if not match:
    print(f"Could not extract time information from filename: {input_file} using regex {file_regex}")
    sys.exit(1)
    
init = dt.datetime.strptime(match.group(1), '%Y%m%d_%H%M')
lead = int(match.group(2))
valid = init + dt.timedelta(hours=lead)

print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

print(valid)
print(init)
print(lead)

attrs = {
        'valid': valid.strftime("%Y%m%d_%H%M%S"),
        'init':  init.strftime("%Y%m%d_%H%M%S"),
        'lead':  str(int(lead)),
        'accum': '00',
        
        'name':      'sept',
        'long_name': 'saturation_equivalent_potential_temperature',
        'level':     'Surface',
        'units':     'K',
        
        'grid': {
            'name': 'Global 0.5 Degree',
            'type' :   'LatLon',
            'lat_ll' : -90.0,
            'lon_ll' : 0.0,
            'delta_lat' :   0.5,
            'delta_lon' :   0.5,
            'Nlat' :      361,
            'Nlon' :      720,
            }
        }

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_ivt_analysis.py

# This script is a combination of two scripts originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR. 
# May 2020
###################################################################################################

import pygrib
import numpy as np
import sys
import os
import re
import datetime as dt
import metpy.calc as mc

###################################################################################################

def ivt(input_file):
    grbs = pygrib.open(input_file)
    g = 9.81    # Setting gravity constant
    print(input_file)
    grbs.rewind()
    
    # Initialize variable arrays 
    levs = [] # Levels 
    q    = [] # Specific humidity
    hgt  = [] # Geopotential height
    temp = [] # Temperature
    u    = [] # u-wind
    v    = [] # v-wind

    # First obtain the levels we will use
    # These are in hPa in the file, so directly compare with user supplied min/max
    levs = sorted(set([grb.level for grb in grbs if float(grb.level) >= float(os.environ.get('IVT_LAYER_MIN_PRESSURE',100.0)) and float(grb.level) <= float(os.environ.get('IVT_LAYER_MAX_PRESSURE',1000.0))]))
    
    # Fill in variable arrays from input file.
    grbs.rewind()
    for grb in grbs:
        if not grb.level in levs:
            continue
        elif np.logical_and('v-' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            v.append(grb.values)        
        elif np.logical_and('u-' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            u.append(grb.values)
        elif np.logical_and('Temperature' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            temp.append(grb.values)
        elif np.logical_and('Geopotential' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            hgt.append(grb.values)
        elif np.logical_and('Specific' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
            q.append(grb.values)
        
    temp = np.array(temp)
    hgt  = np.array(hgt)
    u    = np.array(u)
    v    = np.array(v)
    
    grbs.rewind()
    
    # If we didn't find specific humidity, look for relative humidity.
    if len(q) == 0:
        for grb in grbs:
            if not grb.level in levs:
                continue
            if np.logical_and('Relative' in grb.parameterName,grb.typeOfLevel=='isobaricInhPa'):
                q.append(grb.values)
            
        levs = np.array(levs)
        # Clausius-Clapeyron time
        es = 610.78*np.exp((17.67*(temp-273.15)/(temp-29.65))) # Calculate saturation vapor pressure
        e = es*(np.array(q)/100) # Calculate vapor pressure 
        w = 0.622*es/(levs[:,None,None]*100) # Calculate water vapor
        q = w/(w+1) # Calculate specific humidity
    q = np.array(q)
    
    uv = np.sqrt(u**2+v**2) # Calculate wind
    mflux_total = np.sum(q,axis=0)*(1/g)*np.mean(uv,axis=0)*(np.max(levs)-np.min(levs)) #calculate mass flux
    met_data = mflux_total.copy() #Pass mass flux to be used by MET tools
    print(np.max(met_data))
    #np.save('{}.npy'.format(sys.argv[1]),mflux_total)
    grbs.close()

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = ivt(input_file) #Call function to calculate IVT

met_data = data
met_data = met_data.astype('float64')

# Automatically fill out time information from input file. 
for token in os.path.basename(input_file).replace('-', '_').split('_'):
   if(re.search("[0-9]{8,8}", token)):
       ymd = dt.datetime.strptime(token[0:8],"%Y%m%d")
   elif(re.search("^[0-9]{4}$", token)):
       hh  = int(token[0:2])
   elif(re.search("^[0-9]{3}$", token)):
       day = int(token.replace("", ""))

print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

# GFS Analysis
valid  = ymd  + dt.timedelta(hours=hh)
init = valid

attrs = {
   'valid': valid.strftime("%Y%m%d_%H%M%S"),
   'init':  init.strftime("%Y%m%d_%H%M%S"),
   'lead':  '00',
   'accum': '00',

   'name':      'ivt',
   'long_name': 'integrated_vapor_transport',
   'level':     'Surface',
   'units':     'UNKNOWN',

   'grid': {
       'name': 'Global 0.5 Degree',
       'type' :   'LatLon',
       'lat_ll' : -90.0,
       'lon_ll' : 0.0,
       'delta_lat' :   0.5,
       'delta_lon' :   0.5,
       'Nlat' :      361,
       'Nlon' :      720,
   }
}

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_pv_analysis.py

# This script calculates potential vorticity (PV) from variables found in the GFS analysis model grib2 files. This script is originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR.
# July 2020
###################################################################################################

import sys
import os
import re
import datetime as dt
from metpy import calc as mpcalc
from metpy.units import units
import xarray as xr
import cfgrib

###################################################################################################

def pv(input_file):

    # Vars
    grib_vars = ['t','u','v']

    # Load a list of datasets, one for each variable we want
    ds_list = [cfgrib.open_datasets(input_file,backend_kwargs={'filter_by_keys':{'typeOfLevel':'isobaricInhPa','shortName':v},'indexpath':''}) for v in grib_vars]

    # Flatten the list of lists to a single list of datasets
    ds_flat = [x.sel(isobaricInhPa=x.isobaricInhPa[x.isobaricInhPa>=100.0].values) for ds in ds_list for x in ds]

    # Merge the variables into a single dataset
    ds = xr.merge(ds_flat)

    # Add pressure
    ds['p'] = xr.DataArray(ds.isobaricInhPa.values,dims=['isobaricInhPa'],coords={'isobaricInhPa':ds.isobaricInhPa.values},attrs={'units':'hPa'}).broadcast_like(ds['t'])

    # Calculate potential temperature
    ds['theta'] = mpcalc.potential_temperature(ds['p'].metpy.convert_units('Pa'),ds['t'])

    # Compute baroclinic PV
    ds['pv'] = mpcalc.potential_vorticity_baroclinic(ds['theta'],ds['p'].metpy.convert_units('Pa'),ds['u'],ds['v'],latitude=ds.latitude)/(1.0e-6)

    met_data = ds['pv'].sel(isobaricInhPa=slice(float(os.environ.get('PV_LAYER_MAX_PRESSURE',1000.0)),float(os.environ.get('PV_LAYER_MIN_PRESSURE',100.0)))).mean(axis=0).values

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = pv(input_file) #Call function to calculate PV

met_data = data
met_data = met_data.astype('float64')

print("max", data.max())
print("min", data.min())

# Automatically fill out time information from input file.
for token in os.path.basename(input_file).replace('-', '_').split('_'):
    if(re.search("[0-9]{8,8}", token)):
        ymd = dt.datetime.strptime(token[0:8],"%Y%m%d")
    elif(re.search("^[0-9]{4}$", token)):
        hh  = int(token[0:2])
    elif(re.search("^[0-9]{3}$", token)):
        day = int(token.replace("", ""))


print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

# GFS Analysis
valid  = ymd  + dt.timedelta(hours=hh)
init = valid
#lead, rem = divmod((valid-init).total_seconds(), 3600)


print(valid)
print(init)

attrs = {
    'valid': valid.strftime("%Y%m%d_%H%M%S"),
    'init':  init.strftime("%Y%m%d_%H%M%S"),
    'lead':  '00',
    'accum': '00',
        
    'name':      'pv',
    'long_name': 'potential_vorticity',
    'level':     'Surface',
    'units':     'PV Units',
        
    'grid': {
        'name': 'Global 0.5 Degree',
        'type' :   'LatLon',
        'lat_ll' : -90.0,
        'lon_ll' : 0.0,
        'delta_lat' :   0.5,
        'delta_lon' :   0.5,
        'Nlat' :      361,
        'Nlon' :      720,
    }
}

parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics/gfs_sept_analysis.py

# This script calculates potential vorticity (PV) from variables found in the GFS analysis model grib2 files. This script is originally from Taylor Mandelbaum, SBU. 
# Adjustments have been made by Lindsay Blank, NCAR.
# July 2020
###################################################################################################

import sys
import os
import re
import datetime as dt
from metpy import calc as mpcalc
from metpy.units import units
import xarray as xr
import cfgrib

###################################################################################################

def sept(input_file):

    # Vars
    grib_vars = ['t']

    # Load a list of datasets, one for each variable we want
    ds_list = [cfgrib.open_datasets(input_file,backend_kwargs={'filter_by_keys':{'typeOfLevel':'isobaricInhPa','shortName':v},'indexpath':''}) for v in grib_vars]

    # Flatten the list of lists to a single list of datasets
    ds_flat = [x.sel(isobaricInhPa=x.isobaricInhPa[x.isobaricInhPa>=100.0].values) for ds in ds_list for x in ds]

    # Merge the variables into a single dataset
    ds = xr.merge(ds_flat)

    # Add pressure
    ds['p'] = xr.DataArray(ds.isobaricInhPa.values,dims=['isobaricInhPa'],coords={'isobaricInhPa':ds.isobaricInhPa.values},attrs={'units':'hPa'}).broadcast_like(ds['t'])

    # Calculate saturation equivalent potential temperature
    ds['sept'] = mpcalc.saturation_equivalent_potential_temperature(ds['p'].metpy.convert_units('Pa'),ds['t'])

    met_data = ds['sept'].sel(isobaricInhPa=slice(float(os.environ.get('SEPT_LAYER_MAX_PRESSURE',1000.0)),float(os.environ.get('SEPT_LAYER_MIN_PRESSURE',100.0)))).mean(axis=0).values

    return met_data

###################################################################################################

input_file = os.path.expandvars(sys.argv[1])

data = sept(input_file) #Call function to calculate PV

met_data = data
met_data = met_data.astype('float64')

print("max", data.max())
print("min", data.min())

# Automatically fill out time information from input file.
for token in os.path.basename(input_file).replace('-', '_').split('_'):
    if(re.search("[0-9]{8,8}", token)):
        ymd = dt.datetime.strptime(token[0:8],"%Y%m%d")
    elif(re.search("^[0-9]{4}$", token)):
        hh  = int(token[0:2])
    elif(re.search("^[0-9]{3}$", token)):
        day = int(token.replace("", ""))


print("Data Shape: " + repr(met_data.shape))
print("Data Type:  " + repr(met_data.dtype))

# GFS Analysis
valid  = ymd  + dt.timedelta(hours=hh)
init = valid
#lead, rem = divmod((valid-init).total_seconds(), 3600)


print(valid)
print(init)

attrs = {
    'valid': valid.strftime("%Y%m%d_%H%M%S"),
    'init':  init.strftime("%Y%m%d_%H%M%S"),
    'lead':  '00',
    'accum': '00',
        
    'name':      'sept',
    'long_name': 'saturation_equivalent_potential_temperature',
    'level':     'Surface',
    'units':     'K',
        
    'grid': {
        'name': 'Global 0.5 Degree',
        'type' :   'LatLon',
        'lat_ll' : -90.0,
        'lon_ll' : 0.0,
        'delta_lat' :   0.5,
        'delta_lon' :   0.5,
        'Nlat' :      361,
        'Nlon' :      720,
    }
}

Running METplus

This use case can be run two ways:

1) Passing in TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics.conf, then a user-specific system configuration file:

run_metplus.py \
/path/to/METplus/parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics.conf \
/path/to/user_system.conf
  1. Modifying the configurations in parm/metplus_config, then passing in TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics.conf:

    run_metplus.py \
    /path/to/METplus/parm/use_cases/model_applications/medium_range/TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multiple_Diagnostics.conf
    

The former method is recommended. Whether you add them to a user-specific configuration file or modify the metplus_config files, the following variables must be set correctly:

  • INPUT_BASE - Path to directory where sample data tarballs are unpacked (See Datasets section to obtain tarballs). This is not required to run METplus, but it is required to run the examples in parm/use_cases

  • OUTPUT_BASE - Path where METplus output will be written. This must be in a location where you have write permissions

  • MET_INSTALL_DIR - Path to location where MET is installed locally

and for the [exe] section, you will need to define the location of NON-MET executables. If the executable is in the user’s path, METplus will find it from the name. If the executable is not in the path, specify the full path to the executable here (i.e. CONVERT = /usr/bin/convert) The following executables are required for performing series analysis use cases:

Example User Configuration File:

[dir]
INPUT_BASE = /path/to/sample/input/data
OUTPUT_BASE = /path/to/output/dir
MET_INSTALL_DIR = /path/to/met-X.Y

[exe]
CONVERT = /path/to/convert

Expected Output

A successful run will output the following both to the screen and to the logfile:

INFO: METplus has successfully finished running.

Refer to the value set for OUTPUT_BASE to find where the output data was generated. Output for this use case will be found in subdirectories of the ‘series_analysis_lead’ directory (relative to OUTPUT_BASE):

  • series_animate

  • series_F090

  • series_F096

  • series_F102

  • series_F108

  • series_F114

The series_animate directory contains the animations of the series analysis in .gif format for all variable, level, and statistics combinations:

series_animate_<varname>_<level>_<stat>.gif

The series_FHHH directories contains files that have the following format:

ANLY_FILES_FHHH

FCST_ASCII_FILES_FHHH

series_FHHH_<varname>_<level>_<stat>.png

series_FHHH_<varname>_<level>_<stat>.ps

series_FHHH_<varname>_<level>_<stat>.nc

Where:

HHH is the forecast hour/lead time in hours

varname is the variable of interest, as specified in the METplus series_by_lead_all_fhrs config file

level is the level of interest, as specified in the METplus series_by_lead_all_fhrs config file

stat is the statistic of interest, as specified in the METplus series_by_lead_all_fhrs config file.

Keywords

Note

  • TCPairsToolUseCase

  • SeriesByLeadUseCase

  • TCStatToolUseCase

  • RegridDataPlaneToolUseCase

  • PyEmbedIngestToolUseCase

  • MediumRangeAppUseCase

  • SeriesAnalysisUseCase

  • GRIB2FileUseCase

  • FeatureRelativeUseCase

  • SBUOrgUseCase

  • DiagnosticsUseCase

  • RuntimeFreqUseCase

  • TropicalCycloneUseCase

Navigate to the METplus Quick Search for Use Cases page to discover other similar use cases.

sphinx_gallery_thumbnail_path = ‘_static/medium_range-TCStat_SeriesAnalysis_fcstGFS_obsGFS_FeatureRelative_SeriesByLead_PyEmbed_Multivariate_Diagnostics.png’

Total running time of the script: (0 minutes 0.000 seconds)

Gallery generated by Sphinx-Gallery