UserScript and SeriesAnalysis: Compute Zonal Mean Bias and Create Plots for Temperature and Wind

model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias.py

Scientific Objective

Many common modes of variability in the troposphere have stratospheric teleconnection pathways. Thus, the stratosphere can be a source of predictability for surface weather. However, the predictive skill gained from the stratosphere can be limited by biases in the representation of these stratospheric processes. This use case investigates stratospheric biases in temperature and wind. Model biases are plotted for a month time period based on latitude and pressure between 100 and 1 hPa.

In addition, this use case also demonstrates how to read semi-structured grids into MET. Specifically zonal mean data (on a grid of latitude and pressure) is read into Series-Analysis in the second step of this use case.

Version Added

METplus version 6.0

Datasets

Forecast: GFS 24 hour forecasts for February 2018

Observation: ERA reanalysis for February 2018

Climatology: None

Location: All of the input data required for this use case can be found in a sample data tarball. Each use case category will have one or more sample data tarballs. It is only necessary to download the tarball with the use case’s dataset and not the entire collection of sample data. Click here to access the METplus releases page and download sample data for the appropriate release: https://github.com/dtcenter/METplus/releases This tarball should be unpacked into the directory that you will set the value of INPUT_BASE. See Running METplus section for more information.

METplus Components

This use case calls UserScript first, Series-Analysis, and then UserScript a second time. METcalcpy, METplotpy, and METdataio are needed for this use case to run. The METcalcpy scripts accessed include the following:

  • metcalcpy/pre_processing/directional_means.py

The METplotpy scripts accessed include the following:

  • metplotpy/stratosphere_diagnostics/stratosphere_plots.py

The METdataio scripts accessed include the following:

  • METreadnc/util/read_netcdf.py

METplus Workflow

Beginning time (INIT_BEG): 02-01-2018

End time (INIT_END): 02-28-2018

Increment between beginning and end times (INIT_INCREMENT): 30 days

Sequence of forecast leads to process (LEAD_SEQ): 24 hours

This use case does not loop. The two calls to UserScript are run once. Series- Analysis is also run once. The first call to UserScript runs zonal_mean_driver.py. This script computes zonal and meridional mean data from the directional_means program in METplotpy. Then, Series-Analysis is run on the zonal mean output to compute continuous statistics. Finally, the second call to UserScript runs bias_plot_driver.py which creates bias plots for temperature and wind.

METcalcpy 3.0.0 or higher, METplotpy 3.0.0 or higher, and METdataio 2.1 or higher are needed for this use case.

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, i.e. parm/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias.conf

[config]

# Documentation for this use case can be found at
# https://metplus.readthedocs.io/en/latest/generated/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias.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 = UserScript(obs_means), UserScript(fcst_means), SeriesAnalysis(sa_ut), UserScript(bias_plot)


###
# 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 = VALID

VALID_TIME_FMT = %Y%m%d
VALID_BEG = 20180201
VALID_END = 20180228
VALID_INCREMENT = 30d
LEAD_SEQ = 24

USER_SCRIPT_RUNTIME_FREQ = RUN_ONCE

LOOP_ORDER = processes


###
# User Environment Variables
###
[user_env_vars]
OBS_INPUT_FILE_NAME = {INPUT_BASE}/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/ERA_2018_02.nc
FCST_INPUT_FILE_NAME = {INPUT_BASE}/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/GFS_2018_02_024h.nc

OUTPUT_DIR = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias

PLOT_U_INPUT_FILE = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/SeriesAnalysis/zonal_mean_U_stats_2018_02.nc
PLOT_U_BIAS_VAR = series_cnt_ME
PLOT_U_OBAR_VAR = series_cnt_OBAR
PLOT_U_LEVELS = 0,10,20,30,40,50,60,70,80,90,100
PLOT_U_TITLE = GFS vs ERA 24h Forecast Zonal Mean Wind Bias (ME) 02/2018
PLOT_U_OUTPUT_DIR = {OUTPUT_DIR}/plots
PLOT_U_OUTPUT_FILE = GFS_ERA_ME_2018_02_zonal_mean_U.png

PLOT_T_INPUT_FILE = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/SeriesAnalysis/zonal_mean_T_stats_2018_02.nc
PLOT_T_BIAS_VAR = series_cnt_ME
PLOT_T_OBAR_VAR = series_cnt_OBAR
PLOT_T_LEVELS = 200,210,220,230,240,250,260,270,280,290,300
PLOT_T_TITLE = GFS vs ERA 24h Forecast Zonal Mean Temperature Bias (ME) 02/2018
PLOT_T_OUTPUT_DIR = {OUTPUT_DIR}/plots
PLOT_T_OUTPUT_FILE = GFS_ERA_ME_2018_02_zonal_mean_T.png


###
# Zonal Mean UserScript Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#userscript
###
[obs_means]
USER_SCRIPT_COMMAND = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/zonal_mean_driver.py obs time


[fcst_means]
USER_SCRIPT_COMMAND = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/zonal_mean_driver.py fcst time


###
# Series Analysis Settings
###
[sa_ut]
SERIES_ANALYSIS_RUNTIME_FREQ = RUN_ONCE

FCST_SERIES_ANALYSIS_INPUT_DIR = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/FCST
FCST_SERIES_ANALYSIS_INPUT_TEMPLATE = FCST_zonal_mean_U_T_{valid?fmt=%Y%m}*_{valid?fmt=%H%M%S}.nc
FCST_SERIES_ANALYSIS_INPUT_DATATYPE = PYTHON_NUMPY

OBS_SERIES_ANALYSIS_INPUT_DIR = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/OBS
OBS_SERIES_ANALYSIS_INPUT_TEMPLATE = OBS_zonal_mean_U_T_{valid?fmt=%Y%m}*_{valid?fmt=%H%M%S}.nc
OBS_SERIES_ANALYSIS_INPUT_DATATYPE = PYTHON_NUMPY

SERIES_ANALYSIS_CLIMO_MEAN_INPUT_DIR =
SERIES_ANALYSIS_CLIMO_MEAN_INPUT_TEMPLATE =

SERIES_ANALYSIS_CLIMO_STDEV_INPUT_DIR =
SERIES_ANALYSIS_CLIMO_STDEV_INPUT_TEMPLATE =

SERIES_ANALYSIS_OUTPUT_DIR = {OUTPUT_BASE}/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/SeriesAnalysis

SERIES_ANALYSIS_DESC =

SERIES_ANALYSIS_CAT_THRESH =

SERIES_ANALYSIS_VLD_THRESH =

SERIES_ANALYSIS_BLOCK_SIZE =

SERIES_ANALYSIS_CTS_LIST =

SERIES_ANALYSIS_REGRID_TO_GRID =
SERIES_ANALYSIS_REGRID_METHOD =
SERIES_ANALYSIS_REGRID_WIDTH =
SERIES_ANALYSIS_REGRID_VLD_THRESH =
SERIES_ANALYSIS_REGRID_SHAPE =

SERIES_ANALYSIS_RUN_ONCE_PER_STORM_ID = False

SERIES_ANALYSIS_IS_PAIRED = False

SERIES_ANALYSIS_CONFIG_FILE = {PARM_BASE}/met_config/SeriesAnalysisConfig_wrapped

SERIES_ANALYSIS_OUTPUT_STATS_CNT = TOTAL, ME, RMSE, FBAR, OBAR

MODEL = GFS

OBTYPE = ERA

SERIES_ANALYSIS_OUTPUT_TEMPLATE = zonal_mean_{fcst_level}_stats_2018_02.nc

FCST_SERIES_ANALYSIS_VAR1_NAME = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/read_met_axis_mean.py MET_PYTHON_INPUT_ARG u lat
# This level is used as a label only
FCST_SERIES_ANALYSIS_VAR1_LEVELS = U 

OBS_SERIES_ANALYSIS_VAR1_NAME = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/read_met_axis_mean.py MET_PYTHON_INPUT_ARG u lat


FCST_SERIES_ANALYSIS_VAR2_NAME = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/read_met_axis_mean.py MET_PYTHON_INPUT_ARG T lat
FCST_SERIES_ANALYSIS_VAR2_LEVELS = T

OBS_SERIES_ANALYSIS_VAR2_NAME = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/read_met_axis_mean.py MET_PYTHON_INPUT_ARG T lat


###
# UserScript Bias Plot Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#userscript
###
[bias_plot]
USER_SCRIPT_COMMAND = {PARM_BASE}/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/bias_plot_driver.py

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

SeriesAnalysisConfig_wrapped
////////////////////////////////////////////////////////////////////////////////
//
// 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}
   ${METPLUS_FCST_CLIMO_MEAN_DICT}
   ${METPLUS_FCST_CLIMO_STDEV_DICT}
}
obs = {
   ${METPLUS_OBS_FILE_TYPE}
   ${METPLUS_OBS_CAT_THRESH}
   //field = [
   ${METPLUS_OBS_FIELD}
   ${METPLUS_OBS_CLIMO_MEAN_DICT}
   ${METPLUS_OBS_CLIMO_STDEV_DICT}
}

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

//
// 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}

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

//
// Gradient statistics
// May be set separately in each "obs.field" entry
//
//gradient = {
${METPLUS_GRADIENT_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_TIME_OFFSET_WARNING}
${METPLUS_MET_CONFIG_OVERRIDES}

Python Embedding

This use case uses a Python embedding script to read in the semi-structured zonal mean data to Series-Analysis. Inputs to this script include the filename to be read in, variable name, and the axis over which the mean is taken. The script returns a numpy array containing the zonal mean data (semi-structured grid).

read_met_axis_mean.py
import os
import sys
import numpy as np
import datetime as dt
import netCDF4 as nc4


if len(sys.argv) != 4:
    print("Must specify exactly one input file, variable name, and summary axis (lat, lon, latlon).")
    sys.exit(1)

# Read the input file as the first argument
input_file = os.path.expandvars(sys.argv[1])
var_name   = sys.argv[2]
axis       = sys.argv[3]

# Read the data
f = nc4.Dataset(input_file, 'r')
data = np.float64(f.variables[var_name][:])
met_data = np.transpose(data).copy()

if axis == "lon":
   lats   = list()
   lons   = list(np.float64(f.variables["lon"][:]))
elif axis == "lat":
   lats   = list(np.float64(f.variables["lat"][:]))
   lons   = list()

pres = list(list(np.float64(f.variables["pres"][:])))
times  = list()

lead_ma = f.variables["lead_time"][:]
lead = lead_ma.__int__()
vtime = f.variables["time"]
cur_date = nc4.num2date(vtime[:], vtime.units, vtime.calendar)
init = cur_date - dt.timedelta(hours=lead)
accum     = "00"

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

   'name':      var_name,
   'long_name': str(getattr(f.variables[var_name], "long_name")),
   'level':     axis + "_mean",
   'units':     str(getattr(f.variables[var_name], "units")),

   'grid': {
     'type'   : "SemiLatLon",
     'name'   : axis + "_mean",
     'lats'   : lats,
     'lons'   : lons,
     'levels' : pres,
     'times'  : times
   }
}

print("Attributes: " + repr(attrs))

User Scripting

This use case runs both zonal_mean_driver.py and bias_plot_driver.py. The zonal mean driver takes an input netCDF file and the time variable. Then it computes zonal and meridional means for u and T from directional_means.py in METcalcpy. It writes the zonal mean data to output netCDF files.

The bias plot driver reads output netCDF files from Series-Analysis and creates plots of the bias over latitude and pressure level. Inputs to both of the Python scripts can be found in the [user_env_vars] section of the UserScript_fcstGFS_obsERA_StratosphereBias.conf file

parm/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/zonal_mean_driver.py
#!/usr/bin/env python3

"""
Create meridonial mean statistics

"""
import os
import sys
import logging
import yaml
import xarray as xr  # http://xarray.pydata.org/
import metcalcpy.util.read_env_vars_in_config as readconfig
import metcalcpy.pre_processing.directional_means as directional_means
import METreadnc.util.read_netcdf as read_netcdf


def main():
    """
    Read arguments
    """
    inlabel = sys.argv[1].upper()
    timevar = sys.argv[2]

    """
    Read METplus environment variables
    """
    print('Reading Input Environment Variables')
    input_file = [os.environ[inlabel+'_INPUT_FILE_NAME']]
    output_dir = os.environ['OUTPUT_DIR']
    full_output_dir = os.path.join(output_dir,inlabel)

    """
    Setup logging
    """
    #logfile = "meridonial_mean.log"
    #logging_level = os.environ.get("LOG_LEVEL","logging.INFO")
    #logging.basicConfig(stream=logfile, level=logging_level)

    """
    Read dataset
    """
    print('Reading input data')
    file_reader = read_netcdf.ReadNetCDF()
    ds = file_reader.read_into_xarray(input_file)[0]
    ds = ds.rename({'lat':'latitude',
                    'lon':'longitude'})

    print('Computing Zonal means')
    uzm = directional_means.zonal_mean(ds.u)
    uzm = uzm.assign_coords(lat=("latitude",ds.latitude.values[:,0]))
    Tzm = directional_means.zonal_mean(ds.T)
    Tzm = Tzm.assign_coords(lat=("latitude",ds.latitude.values[:,0]))

    """
    Write output files if desired, by first creating a directory
    """
    print('Writing output zonal mean files')
    if not os.path.exists(full_output_dir):
        os.makedirs(full_output_dir)

    datetimeindex = ds.indexes[timevar].to_datetimeindex()
    out_ds = uzm.to_dataset()
    out_ds.u.attrs = ds.u.attrs
    out_ds['T'] = Tzm
    out_ds.T.attrs = ds.T.attrs
    for i in range(len(datetimeindex)):
        cur_date = datetimeindex[i]
        output_file = os.path.join(full_output_dir,inlabel+'_zonal_mean_U_T_'+cur_date.strftime('%Y%m%d_%H%M%S')+'.nc')
        out_write = out_ds.isel(time=i)
        # Add lead time as a variable
        if inlabel == 'OBS':
            out_write = out_write.assign(lead_time=[0.0])
        elif inlabel == 'FCST':
            #Grab Forecast Lead
            out_write = out_write.assign(lead_time=ds.lead_time[i])
        out_write.to_netcdf(output_file, 'w')



if __name__ == '__main__':
    main()
parm/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias/bias_plot_driver.py
#!/usr/bin/env python3

"""
Creates a bias plot
"""
import os
import METreadnc.util.read_netcdf as read_netcdf
from metplotpy.contributed.stratosphere_diagnostics.stratosphere_plots import plot_zonal_bias


def main():
    """
    Read METplus environment variables
    """
    print('Reading Input Environment Variables')
    infile_u = [os.environ['PLOT_U_INPUT_FILE']]
    invar_u = os.environ['PLOT_U_BIAS_VAR']
    omvar_u = os.environ['PLOT_U_OBAR_VAR']
    plot_levels_u_str = os.environ['PLOT_U_LEVELS'].split(',')
    plot_levels_u = [int(pp) for pp in plot_levels_u_str]
    plot_title_u = os.environ['PLOT_U_TITLE']
    output_dir_u = os.environ['PLOT_U_OUTPUT_DIR']
    output_file_u = os.environ.get('PLOT_U_OUTPUT_FILE','bias_plot.png')
    plot_output_file_u = os.path.join(output_dir_u,output_file_u)

    infile_t = [os.environ['PLOT_T_INPUT_FILE']]
    invar_t = os.environ['PLOT_T_BIAS_VAR']
    omvar_t = os.environ['PLOT_T_OBAR_VAR']
    plot_levels_t_str = os.environ['PLOT_T_LEVELS'].split(',')
    plot_levels_t = [int(pp) for pp in plot_levels_t_str]
    plot_title_t = os.environ['PLOT_T_TITLE']
    output_dir_t = os.environ['PLOT_T_OUTPUT_DIR']
    output_file_t = os.environ.get('PLOT_T_OUTPUT_FILE','bias_plot.png')
    plot_output_file_t = os.path.join(output_dir_t,output_file_t)

    """
    Make plot output directory if it doesn't exist
    """
    if not os.path.exists(output_dir_u):
        os.makedirs(output_dir_u)

    if not os.path.exists(output_dir_t):
        os.makedirs(output_dir_t)

    """
    Read dataset
    """
    print('Reading Zonal Mean U Bias File: '+infile_u[0])
    file_reader_u = read_netcdf.ReadNetCDF()
    dsu = file_reader_u.read_into_xarray(infile_u)[0]
    bias_u = dsu[invar_u].values
    lats_u = dsu['lat'].values
    obar_u = dsu[omvar_u].values
    levels_u = dsu['level'].values

    print('Reading Zonal Mean T Bias File: '+infile_t[0])
    file_reader_t = read_netcdf.ReadNetCDF()
    dst = file_reader_t.read_into_xarray(infile_t)[0]
    bias_t = dst[invar_t].values
    lats_t = dst['lat'].values
    obar_t = dst[omvar_t].values
    levels_t = dst['level'].values

    """
    Create Bias Plots
    """
    print('Plotting Zonal Mean U Bias')
    plot_zonal_bias(lats_u,levels_u,bias_u,obar_u,plot_output_file_u,plot_title_u,plot_levels_u)
    print('Plotting Zonal Mean T Bias')
    plot_zonal_bias(lats_t,levels_t,bias_t,obar_t,plot_output_file_t,plot_title_t,plot_levels_t)


if __name__ == '__main__':
    main()

Running METplus

Pass the use case configuration file to the run_metplus.py script along with any user-specific system configuration files if desired:

run_metplus.py /path/to/METplus/parm/use_cases/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratosphereBias.conf /path/to/user_system.conf

See Running METplus for more information.

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 {OUTPUT_BASE}/model_applications/s2s_stratosphere/UserScript_fcstGFS_obsERA_StratospherePolar The output includes the netCDF zonal mean files for the forecast and observations, netCDF output from Series-Analysis, and bias plots. There are two bias plots output:

* plots/GFS_ERA_ME_2018_02_zonal_mean_T.png
* plots/GFS_ERA_ME_2018_02_zonal_mean_U.png

The statistics output from Series-Analysis output is two netCDF files, one for temperature and one for wind:

* SeriesAnalysis/zonal_mean_T_stats_2018_02.nc
* SeriesAnalysis/zonal_mean_U_stats_2018_02.nc

There are 7 variable fields present in the Series-Analysis output netCDF files (not including the lat/lon fields). Those variables are:

* level(level)
* n_series
* series_cnt_TOTAL(lat, level)
* series_cnt_ME(lat, level)
* series_cnt_RMSE(lat, level)
* series_cnt_FBAR(lat, level)
* series_cnt_OBAR(lat, level)

Text files with a listing of the files input to Series-Analysis are also output:

* SeriesAnalysis/series_analysis_files_fcst_init_ALL_valid_ALL_lead_ALL.txt
* SeriesAnalysis/series_analysis_files_obs_init_ALL_valid_ALL_lead_ALL.txt

The zonal mean output includes 28 files for the forecast and observations, one for each day. The file format for February 1 is:

* FCST/FCST_zonal_mean_U_T_20180201_000000.nc
* OBS/OBS_zonal_mean_U_T_20180201_000000.nc

There are 4 variable fields present in the zonal mean netCDF file (not including the latitude and pressure fields). Those variables are:

* time
* u(pres, latitude)
* T(pres, latitude)
* lead_time

Keywords

Note

  • UserScriptUseCase

  • S2SAppUseCase

  • S2SStratosphereAppUseCase

  • UserScriptUseCase

  • SeriesAnalysisUseCase

  • METdataioUseCase

  • METcalcpyUseCase

  • METplotpyUseCase

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

sphinx_gallery_thumbnail_path = ‘_static/s2s_stratosphere-UserScript_fcstGFS_obsERA_StratosphereBias.png’

Gallery generated by Sphinx-Gallery