UserScript: Calculate Blocking for the ERA

model_applications/s2s_mid_lat/UserScript_obsERA_obsOnly_Blocking.py

Scientific Objective

Atmospheric blocking is associated with extreme weather events. This use case computes atmospheric blocking events using the methodology in Miller & Wang (2019, 2022), which identifies blocks from 500 mb height. Various studies (Masato et al. 2013; Kitano and Yamada 2016) have suggested that 500 mb height produces a similar blocking climatology as results from potential temperature on a 2-PVU surface.

The methodology in Miller & Wang (2019, 2022) first computes the Central Blocking Latitude (CBL) or storm track. Allowing for an offset north and south of the storm track, reversals in geopotential height are then identified as Instantaneously Blocked Longitudes (IBLs). These IBLs are grouped when consective longitudes are blocked (GIBLs) and then blocks are identified by applying thresholds to ensure the large-scale, quasi-stationary characteristics of blocking anticyclones are met.

This use case is a simplified version of the UserScript_fcstGFS_obs_ERA_Blocking use case. While that use case evaluates blocking for both the model and observations, this case shows how to run the blocking calculation on observations only. The setup is simpler and requires fewer MET runs than the UserScript_fcstGFS_obs_ERA_Blocking use case. The original code for computing blocking came from Douglas Miller.

  • Miller, D. E., and Z. Wang, 2019a: Skillful seasonal prediction of Eurasian winter blocking and extreme temperature frequency. Geophys. Res. Lett., 46, 11 530–11 538, https://doi.org/10.1029/2019GL085035.

  • Miller, D. E., and Z. Wang, 2022: Northern Hemisphere Winter Blocking: Differing Onset Mechanisms across regions. J. Atmos. Sci., 79, 1291-1309, https://doi.org/10.1175/JAS-D-21-0104.1.

  • Masato, G., B. J. Hoskins, and T. J. Woollings, 2013: Winter and summer Northern Hemisphere blocking in CMIP5 models. J. Climate, 26, 7044–7059, https://doi.org/10.1175/JCLI-D-12-00466.1.

  • Kitano, Y., and T. J. Yamada, 2016: Relationship between atmospheric blocking and cold day extremes in current and RCP8.5 future climate conditions over Japan and the surrounding area. Atmos. Sci. Lett., 17, 616–622, https://doi.org/10.1002/asl.711.

Version Added

METplus version 4.0.0

Datasets

Forecast: None

Observation: ERA Reanlaysis 500 mb height for DJF 1979 - 2017

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 once. There are 4 optional pre-processing steps, RegridDataPlane and 3 calls to PcpCombine. METcalcpy and METplotpy are also required to run this use case. The METcalcpy scripts accessed include the following:

  • metcalcpy/contributed/blocking_weather_regime/Blocking.py

  • metcalcpy/contributed/blocking_weather_regime/Blocking_WeatherRegime_util.py

The METplopty scrips accessed include the following:

  • metplotpy/contributed/blocking_s2s/CBL_plot.py

  • metplotpy/contributed/blocking_s2s/plot_blocking.py

METplus Workflow

Beginning time (VALID_BEG): 12-01-1979

End time (VALID_END): 02-28-2017

Increment between beginning and end times (VALID_INCREMENT): 1 day

Sequence of forecast leads to process (LEAD_SEQ): 0 hour

This use case does not loop. It runs UserScript once over the entire time period. The UserScript runs the blocking driver script which performs the blocking calculation. This calculation is divided up into steps, which the user can select by setting STEPS_OBS in the [user_env_vars] section of the configuration. More information on the steps and how the calculation proceeds is given in the User Scripting section below.

The four optional pre-processings steps loop by valid time with different timing settings needed used for the different steps. These steps are turned off due to data size and processing time. The first optional step calls Regrid-Data-Plane to regrid the data to a 1 degree latitude/longitude grid. The second calls PCP-Combine to compute daily means of 500 mb height. The third calls PCP-Combine to compute a 5 day running mean, while the last uses PCP-Combine to compute anomalies by subtracting the 5 day running mean from the daily mean. These omitted steps can be turned on by using the PROCESS_LIST that is commented out:

PROCESS_LIST = RegridDataPlane(regrid_obs), PcpCombine(daily_mean_obs), PcpCombine(running_mean_obs), PcpCombine(anomaly_obs), UserScript(script_blocking)

Settings for the optional pre-processing steps can be found in the respective sections of the configuration, regrid_obs, daily_mean_obs, running_mean_obs, and anomaly_obs. Data is not provided in the tarball to run these steps, but the configurations are provided as an example of how to set up these steps.

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_mid_lat/UserScript_obsERA_obsOnly_Blocking.conf

[config]

# Documentation for this use case can be found at
# https://metplus.readthedocs.io/en/latest/generated/model_applications/s2s_mid_lat/UserScript_obsERA_obsOnly_Blocking.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
###

# All steps, including pre-processing:
# PROCESS_LIST = RegridDataPlane(regrid_obs), PcpCombine(daily_mean_obs), PcpCombine(running_mean_obs), PcpCombine(anomaly_obs), UserScript(script_blocking)
# Only Blocking Analysis script for the observations

PROCESS_LIST = UserScript(script_blocking)


###
# 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%H
VALID_BEG = 1979120100
VALID_END = 2017022800
VALID_INCREMENT = 86400

LEAD_SEQ = 0

# Only Process DJF
SKIP_TIMES = "%m:begin_end_incr(3,11,1)", "%m%d:0229"

# Run the obs data
# A variable set to be used in the pre-processing steps
OBS_RUN = True


###
# RegridDataPlane(regrid_obs) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#regriddataplane
###

# Regrid the observations to 1 degree using regrid_data_plane
[regrid_obs]

VALID_END = 2017022818
VALID_INCREMENT = 21600

OBS_REGRID_DATA_PLANE_RUN = {OBS_RUN}

OBS_DATA_PLANE_ONCE_PER_FIELD = False

OBS_REGRID_DATA_PLANE_VAR1_INPUT_FIELD_NAME = Z
OBS_REGRID_DATA_PLANE_VAR1_INPUT_LEVEL = P500
OBS_REGRID_DATA_PLANE_VAR1_OUTPUT_FIELD_NAME = Z500

REGRID_DATA_PLANE_VERIF_GRID = latlon 360 90 89 0 -1.0 1.0

REGRID_DATA_PLANE_METHOD = BILIN

REGRID_DATA_PLANE_WIDTH = 2

OBS_REGRID_DATA_PLANE_INPUT_DIR = /gpfs/fs1/collections/rda/data/ds627.0/ei.oper.an.pl
OBS_REGRID_DATA_PLANE_OUTPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Regrid

OBS_REGRID_DATA_PLANE_INPUT_TEMPLATE = {valid?fmt=%Y%m}/ei.oper.an.pl.regn128sc.{valid?fmt=%Y%m%d%H}
OBS_REGRID_DATA_PLANE_OUTPUT_TEMPLATE = {valid?fmt=%Y%m%d}/Z500_6hourly_{init?fmt=%Y%m%d%H}_NH.nc


###
# PCPCombine(daily_mean_obs) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#pcpcombine
###

# Perform a sum over the 4 daily times that have been regridded using pcp_combine
# 00, 06, 12, 18 UTC
[daily_mean_obs]

VALID_BEG = 1979120118
VALID_END = 2017022818

OBS_PCP_COMBINE_RUN = {OBS_RUN}

OBS_PCP_COMBINE_METHOD = DERIVE
OBS_PCP_COMBINE_STAT_LIST = MEAN

OBS_PCP_COMBINE_INPUT_ACCUMS = 6
OBS_PCP_COMBINE_INPUT_NAMES = Z500
OBS_PCP_COMBINE_INPUT_LEVELS = "(*,*)"
OBS_PCP_COMBINE_INPUT_OPTIONS = convert(x) = x / 9.81; set_attr_valid = "{valid?fmt=%Y%m%d_%H%M%S?shift=-64800}";

OBS_PCP_COMBINE_OUTPUT_ACCUM = 24
OBS_PCP_COMBINE_DERIVE_LOOKBACK = 24

OBS_PCP_COMBINE_OUTPUT_NAME = Z500

OBS_PCP_COMBINE_INPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Regrid
OBS_PCP_COMBINE_OUTPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Daily

OBS_PCP_COMBINE_INPUT_TEMPLATE = {valid?fmt=%Y%m%d}/Z500_6hourly_{valid?fmt=%Y%m%d%H}_NH.nc
OBS_PCP_COMBINE_OUTPUT_TEMPLATE = Z500_daily_{valid?fmt=%Y%m%d?shift=-64800}_NH.nc


###
# PCPCombine(running_mean_obs) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#pcpcombine
###

# Perform a 5 day running mean on the data using pcp_combine

[running_mean_obs]

# Add the first/last 2 days to the skip times to compute the running mean
SKIP_TIMES = "%m:begin_end_incr(3,11,1)", "%m%d:1201,1202,1203,1204,0229"

OBS_PCP_COMBINE_RUN = {OBS_RUN}

OBS_PCP_COMBINE_METHOD = DERIVE
OBS_PCP_COMBINE_STAT_LIST = MEAN

OBS_PCP_COMBINE_INPUT_ACCUMS = 24
OBS_PCP_COMBINE_INPUT_NAMES = Z500
OBS_PCP_COMBINE_INPUT_LEVELS = "(*,*)"
OBS_PCP_COMBINE_INPUT_OPTIONS = set_attr_valid = "{valid?fmt=%Y%m%d_%H%M%S?shift=-172800}";

#  Running mean is 5 days
OBS_PCP_COMBINE_OUTPUT_ACCUM = 120
OBS_PCP_COMBINE_DERIVE_LOOKBACK = 120

OBS_PCP_COMBINE_OUTPUT_NAME = Z500

OBS_PCP_COMBINE_INPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Daily
OBS_PCP_COMBINE_OUTPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Rmean5d

OBS_PCP_COMBINE_INPUT_TEMPLATE = Z500_daily_{valid?fmt=%Y%m%d}_NH.nc
OBS_PCP_COMBINE_OUTPUT_TEMPLATE = Z500_5daymean_{valid?fmt=%Y%m%d?shift=-172800}_NH.nc


###
# PCPCombine(anomaly_obs) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#pcpcombine
###

# Compute anomalies using the daily means and 5 day running mean using pcp_combine

[anomaly_obs]

# Add the first/last 2 days to the skip times to compute the running mean 
SKIP_TIMES = "%m:begin_end_incr(3,11,1)", "%m%d:1201,1202,0227,0228,0229"

OBS_PCP_COMBINE_RUN = {OBS_RUN}

OBS_PCP_COMBINE_METHOD = USER_DEFINED

OBS_PCP_COMBINE_COMMAND = -subtract {OBS_PCP_COMBINE_INPUT_DIR}/Daily/Z500_daily_{valid?fmt=%Y%m%d}_NH.nc {OBS_PCP_COMBINE_INPUT_DIR}/Rmean5d/Z500_5daymean_{valid?fmt=%Y%m%d}_NH.nc -field 'name="Z500"; level="(*,*)";'

OBS_PCP_COMBINE_INPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA
OBS_PCP_COMBINE_OUTPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Anomaly

OBS_PCP_COMBINE_INPUT_TEMPLATE = Z500_daily_{valid?fmt=%Y%m%d}_NH.nc
OBS_PCP_COMBINE_OUTPUT_TEMPLATE = Z500_anomaly_{valid?fmt=%Y%m%d}_NH.nc


# Variables set for the Blocking Analysis
[user_env_vars]
# Steps to Run
OBS_STEPS = CBL+PLOTCBL+IBL+PLOTIBL+GIBL+CALCBLOCKS+PLOTBLOCKS

# Number of Seasons and Days per season that should be available
# The code will fill missing data, but requires the same number of days per
# season for each year.  You may need to omit leap days if February is part of
# the processing
CBL_NUM_SEASONS = 38
IBL_NUM_SEASONS = 38
DAYS_PER_SEASON = 86

# Make the OUTPUT_BASE available to the UserScript
SCRIPT_OUTPUT_BASE = {OUTPUT_BASE}

# Variable Name for the Z500 anomaly data to read in to the blocking python code
OBS_BLOCKING_ANOMALY_VAR = Z500_ANA

# Variable for the Z500 data
OBS_BLOCKING_VAR = Z500

# Number of model grid points used for a moving average
# Must be odd
OBS_SMOOTHING_PTS = 9

# Lat Delta, to allow for offset from the Central Blocking Latitude
OBS_LAT_DELTA = -5,0,5

# Meridional Extent of blocks (NORTH_SOUTH_LIMITS/2)
OBS_NORTH_SOUTH_LIMITS = 30

# Maximum number of grid points between IBLs for everything in between to be included as an IBL
OBS_IBL_DIST = 7

# Number of grid points in and IBL to make a GIBL
OBS_IBL_IN_GIBL = 15

# Number of grid points that must overlap across days for a GIBL
OBS_GIBL_OVERLAP = 10

# Time duration in days needed for a block
OBS_BLOCK_TIME = 5

# Number of grid points a block must travel to terminate
OBS_BLOCK_TRAVEL = 45

# Method to compute blocking.  Currently, the only option is 'PH' for the
# Pelly-Hoskins Method
OBS_BLOCK_METHOD = PH

# Plot Output Directory
BLOCKING_PLOT_OUTPUT_DIR = {OUTPUT_BASE}/s2s_mid_lat/UserScript_obsERA_obsOnly_Blocking/plots/

#CBL plot title and name
OBS_CBL_PLOT_MTHSTR = DJF
OBS_CBL_PLOT_OUTPUT_NAME = ERA_CBL_avg

# IBL plot title and name
OBS_IBL_PLOT_TITLE = DJF ERA Instantaneous Blocked Longitude
OBS_IBL_PLOT_OUTPUT_NAME = ERA_IBL_Freq_DJF

# Blocking plot title and name
OBS_BLOCKING_PLOT_TITLE = DJF ERA Blocking Frequency
OBS_BLOCKING_PLOT_OUTPUT_NAME = ERA_Block_Freq_DJF


###
# UserScript(script_blocking) Settings
# https://metplus.readthedocs.io/en/latest/Users_Guide/wrappers.html#userscript
###

# Run the Blocking Analysis Script
[script_blocking]

# Skip the days on the edges that are not available due to the running mean
SKIP_TIMES = "%m:begin_end_incr(3,11,1)", "%m%d:1201,1202,0227,0228,0229"

# Run the user script once per lead
USER_SCRIPT_RUNTIME_FREQ = RUN_ONCE_PER_LEAD

# Template of filenames to input to the user-script
USER_SCRIPT_INPUT_TEMPLATE = {INPUT_BASE}/model_applications/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Anomaly/Z500_anomaly_{valid?fmt=%Y%m%d}_NH.nc,{INPUT_BASE}/model_applications/s2s_mid_lat/UserScript_fcstGFS_obsERA_Blocking/ERA/Daily/Z500_daily_{valid?fmt=%Y%m%d}_NH.nc

# Name of the file containing the listing of input files
# The options are OBS_CBL_INPUT, FCST_CBL_INPUT, OBS_IBL_INPUT, and FCST_IBL_INPUT
# *** Make sure the order is the same as the order of templates listed in USER_SCRIPT_INPUT_TEMPLATE
USER_SCRIPT_INPUT_TEMPLATE_LABELS = OBS_CBL_INPUT,OBS_IBL_INPUT

# Command to run the user script with input configuration file
USER_SCRIPT_COMMAND = {METPLUS_BASE}/parm/use_cases/model_applications/s2s_mid_lat/common/Blocking_driver.py

MET Configuration

This case does not use MET configuration files.

Python Embedding

This use case does not use python embedding.

User Scripting

This use case runs the blocking driver. The blocking driver runs the user selected steps of the blocking calculation. These steps are specified in OBS_STEPS in the [user_env_vars] section fo the configuration file in the following format:

OBS_STEPS = CBL+PLOTCBL+IBL+PLOTIBL+GILB+CALCBLOCKS+PLOTBLOCKS

The possible steps are computing the CBLs or Central Blocking Latitude (CBL), plotting CBLs (PLOTCBL), computing Instantaneously Blocked Longitudes (IBL), plotting IBL frequency (PLOTIBL), computing grouped instantaneously blocked longitudes (GIBL), computing blocks (CALCBLOCKS), and plotting the blocking frequency (PLOTBLOCKS). This use case runs all steps although not all of them are required to be run. The CBL, IBL, GIBL, and CALCBLOCKS steps must be run in order as the IBL step requires previously computed CBLs, and GIBLs requires previously computed IBLs. Plotting also requires the associated step to be run (e.g. PLOTCBL requires CBL to be run first). The methodology used in these calculations is described in Miller & Wang (2019, 2022) listed in the Scientific Objective section.

There are many input variables that can be changed for the driver script and blocking calculation. These can be changed and are described in the [user_env_vars] section of the configuration file.

parm/use_cases/model_applications/s2s_mid_lat/common/Blocking_driver.py
#!/usr/bin/env python3

import os
import warnings

import numpy as np

from metcalcpy.contributed.blocking_weather_regime.Blocking import BlockingCalculation
from metcalcpy.contributed.blocking_weather_regime.Blocking_WeatherRegime_util import parse_steps, write_mpr_file
from metplotpy.contributed.blocking_s2s import plot_blocking as pb
from metplotpy.contributed.blocking_s2s.CBL_plot import create_cbl_plot


def main():

    steps_list_fcst,steps_list_obs = parse_steps()

    if not steps_list_obs and not steps_list_fcst:
        warnings.warn('No processing steps requested for either the model or observations,')
        warnings.warn(' nothing will be run')
        warnings.warn('Set FCST_STEPS and/or OBS_STEPS in the [user_env_vars] section to process data')


    ######################################################################
    # Blocking Calculation and Plotting
    ######################################################################
    # Set up the data
    steps_fcst = BlockingCalculation('FCST')
    steps_obs = BlockingCalculation('OBS')

    # Check to see if there is a plot directory
    oplot_dir = os.environ.get('BLOCKING_PLOT_OUTPUT_DIR','')
    if not oplot_dir:
        obase = os.environ['SCRIPT_OUTPUT_BASE']
        oplot_dir = os.path.join(obase,'plots')
    if not os.path.exists(oplot_dir):
        os.makedirs(oplot_dir)

    # Check to see if there is a mpr output directory
    mpr_dir = os.environ.get('BLOCKING_MPR_OUTPUT_DIR','')
    if not mpr_dir:
        obase = os.environ['SCRIPT_OUTPUT_BASE']
        mpr_dir = os.path.join(obase,'mpr')

    # Check to see if CBL's are used from an obs climatology
    use_cbl_obs = os.environ.get('USE_CBL_OBS','False').lower()

    # Get the days per season
    dseasons = int(os.environ['DAYS_PER_SEASON'])

    # Grab the Anomaly (CBL) text files
    obs_cbl_filetxt = os.environ.get('METPLUS_FILELIST_OBS_CBL_INPUT','')
    fcst_cbl_filetxt = os.environ.get('METPLUS_FILELIST_FCST_CBL_INPUT','')

    # Grab the Daily (IBL) text files
    obs_ibl_filetxt = os.environ.get('METPLUS_FILELIST_OBS_IBL_INPUT','')
    fcst_ibl_filetxt = os.environ.get('METPLUS_FILELIST_FCST_IBL_INPUT','')


    # Calculate Central Blocking Latitude
    if "CBL" in steps_list_obs:
        print('Computing Obs CBLs')
        # Read in the list of CBL files
        cbl_nseasons = int(os.environ['CBL_NUM_SEASONS'])
        with open(obs_cbl_filetxt) as ocl:
            obs_infiles = ocl.read().splitlines()
        if obs_infiles[0] == 'file_list':
            obs_infiles = obs_infiles[1:]
        if len(obs_infiles) != (cbl_nseasons*dseasons):
            raise ValueError('Invalid Obs data; each year must contain the same date range to calculate seasonal averages.')
        cbls_obs,lats_obs,lons_obs,mhweight_obs,_ = steps_obs.run_CBL(obs_infiles,cbl_nseasons,dseasons)

    if "CBL" in steps_list_fcst and use_cbl_obs == 'false':
        # Add in step to use obs for CBLS
        print('Computing Forecast CBLs')
        cbl_nseasons = int(os.environ['CBL_NUM_SEASONS'])
        with open(fcst_cbl_filetxt) as fcl:
            fcst_infiles = fcl.read().splitlines()
        if fcst_infiles[0] == 'file_list':
            fcst_infiles = fcst_infiles[1:]
        if len(fcst_infiles) != (cbl_nseasons*dseasons):
            raise ValueError('Invalid Fcst data; each year must contain the same date range to calculate seasonal averages.')
        cbls_fcst,lats_fcst,lons_fcst,mhweight_fcst,_ = steps_fcst.run_CBL(fcst_infiles,cbl_nseasons,dseasons)
    elif ("CBL" in steps_list_fcst) and (use_cbl_obs == 'true'):
        if "CBL" not in steps_list_obs:
            raise ValueError('Must run observed CBLs before using them as a forecast.')
        cbls_fcst = cbls_obs
        lats_fcst = lats_obs
        lons_fcst = lons_obs
        mhweight_fcst = mhweight_obs

    #Plot Central Blocking Latitude
    if "PLOTCBL" in steps_list_obs:
        if "CBL" not in steps_list_obs:
            raise ValueError('Must run observed CBLs before plotting them.')
        print('Plotting Obs CBLs')
        cbl_plot_mthstr = os.environ['OBS_CBL_PLOT_MTHSTR']
        cbl_plot_outname = os.path.join(oplot_dir,os.environ.get('OBS_CBL_PLOT_OUTPUT_NAME','obs_cbl_avg'))
        create_cbl_plot(lons_obs, lats_obs, cbls_obs, mhweight_obs, cbl_plot_mthstr, cbl_plot_outname, 
            do_averaging=True)
    if "PLOTCBL" in steps_list_fcst:
        if "CBL" not in steps_list_fcst:
            raise ValueError('Must run forecast CBLs before plotting them.')
        print('Plotting Forecast CBLs')
        cbl_plot_mthstr = os.environ['FCST_CBL_PLOT_MTHSTR']
        cbl_plot_outname = os.path.join(oplot_dir,os.environ.get('FCST_CBL_PLOT_OUTPUT_NAME','fcst_cbl_avg'))
        create_cbl_plot(lons_fcst, lats_fcst, cbls_fcst, mhweight_fcst, cbl_plot_mthstr, cbl_plot_outname, 
            do_averaging=True)


    # Run IBL
    if "IBL" in steps_list_obs:
        if "CBL" not in steps_list_obs:
            raise ValueError('Must run observed CBLs before running IBLs.')
        print('Computing Obs IBLs')
        ibl_nseasons = int(os.environ['IBL_NUM_SEASONS'])
        with open(obs_ibl_filetxt) as oil:
            obs_infiles = oil.read().splitlines()
        if obs_infiles[0] == 'file_list':
            obs_infiles = obs_infiles[1:]
        if len(obs_infiles) != (ibl_nseasons*dseasons):
            raise ValueError('Invalid Obs data; each year must contain the same date range to calculate seasonal averages.')
        ibls_obs,ibl_time_obs = steps_obs.run_Calc_IBL(cbls_obs,obs_infiles,ibl_nseasons,dseasons)
        daynum_obs = np.arange(0,len(ibls_obs[0,:,0]),1)
    if "IBL" in steps_list_fcst:
        if "CBL" not in steps_list_fcst:
            raise ValueError('Must run forecast CBLs or use observed CBLs before running IBLs.')
        print('Computing Forecast IBLs')
        ibl_nseasons = int(os.environ['IBL_NUM_SEASONS'])
        with open(fcst_ibl_filetxt) as fil:
            fcst_infiles = fil.read().splitlines()
        if fcst_infiles[0] == 'file_list':
            fcst_infiles = fcst_infiles[1:]
        if len(fcst_infiles) != (ibl_nseasons*dseasons):
            raise ValueError('Invalid Fcst data; each year must contain the same date range to calculate seasonal averages.')
        ibls_fcst,ibl_time_fcst = steps_fcst.run_Calc_IBL(cbls_fcst,fcst_infiles,ibl_nseasons,dseasons)
        daynum_fcst = np.arange(0,len(ibls_fcst[0,:,0]),1)

    if "IBL" in steps_list_obs and "IBL" in steps_list_fcst:
        # Print IBLs to output matched pair file
        i_mpr_outdir = os.path.join(mpr_dir,'IBL')
        if not os.path.exists(i_mpr_outdir):
            os.makedirs(i_mpr_outdir)
        modname = os.environ.get('MODEL_NAME','GFS')
        maskname = os.environ.get('MASK_NAME','FULL')
        ibl_outfile_prefix = os.path.join(i_mpr_outdir,'IBL_stat_'+modname)
        cbls_avg = np.nanmean(cbls_obs,axis=0)
        write_mpr_file(ibls_obs,ibls_fcst,cbls_avg,lons_obs,ibl_time_obs,ibl_time_fcst,modname,
            'NA','IBLs','block','Z500','IBLs','block','Z500',maskname,'500',ibl_outfile_prefix)

    # Plot IBLS
    if "PLOTIBL" in steps_list_obs and "PLOTIBL" not in steps_list_fcst:
        if "IBL" not in steps_list_obs:
            raise ValueError('Must run observed IBLs before plotting them.')
        print('Plotting Obs IBLs')
        ibl_plot_title = os.environ.get('OBS_IBL_PLOT_TITLE','Instantaneous Blocked Longitude')
        ibl_plot_outname = os.path.join(oplot_dir,os.environ.get('OBS_IBL_PLOT_OUTPUT_NAME','obs_IBL_Freq'))
        ibl_plot_label1 = os.environ.get('IBL_PLOT_OBS_LABEL','')
        pb.plot_ibls(ibls_obs,lons_obs,ibl_plot_title,ibl_plot_outname,label1=ibl_plot_label1)
    elif "PLOTIBL" in steps_list_fcst and "PLOTIBL" not in steps_list_obs:
        if "IBL" not in steps_list_fcst:
            raise ValueError('Must run forecast IBLs before plotting them.')
        print('Plotting Forecast IBLs')
        ibl_plot_title = os.environ.get('FCST_IBL_PLOT_TITLE','Instantaneous Blocked Longitude')
        ibl_plot_outname = os.path.join(oplot_dir,os.environ.get('FCST_IBL_PLOT_OUTPUT_NAME','fcst_IBL_Freq'))
        ibl_plot_label1 = os.environ.get('IBL_PLOT_FCST_LABEL','')
        pb.plot_ibls(ibls_fcst,lons_fcst,ibl_plot_title,ibl_plot_outname,label1=ibl_plot_label1)
    elif ("PLOTIBL" in steps_list_obs) and ("PLOTIBL" in steps_list_fcst):
        if "IBL" not in steps_list_obs and "IBL" not in steps_list_fcst:
            raise ValueError('Must run forecast and observed IBLs before plotting them.')
        print('Plotting Obs and Forecast IBLs')
        ibl_plot_title = os.environ['IBL_PLOT_TITLE']
        ibl_plot_outname = os.path.join(oplot_dir,os.environ.get('IBL_PLOT_OUTPUT_NAME','IBL_Freq'))
        #Check to see if there are plot legend labels
        ibl_plot_label1 = os.environ.get('IBL_PLOT_OBS_LABEL','Observation')
        ibl_plot_label2 = os.environ.get('IBL_PLOT_FCST_LABEL','Forecast')
        pb.plot_ibls(ibls_obs,lons_obs,ibl_plot_title,ibl_plot_outname,data2=ibls_fcst,lon2=lons_fcst,
            label1=ibl_plot_label1,label2=ibl_plot_label2)


    # Run GIBL
    if "GIBL" in steps_list_obs:
        if "IBL" not in steps_list_obs:
            raise ValueError('Must run observed IBLs before running GIBLs.')
        print('Computing Obs GIBLs')
        gibls_obs = steps_obs.run_Calc_GIBL(ibls_obs,lons_obs)

    if "GIBL" in steps_list_fcst:
        if "IBL" not in steps_list_fcst:
            raise ValueError('Must run Forecast IBLs before running GIBLs.')
        print('Computing Forecast GIBLs')
        gibls_fcst = steps_fcst.run_Calc_GIBL(ibls_fcst,lons_fcst)


    # Calc Blocks
    if "CALCBLOCKS" in steps_list_obs:
        if "GIBL" not in steps_list_obs:
            raise ValueError('Must run observed GIBLs before calculating blocks.')
        print('Computing Obs Blocks')
        block_freq_obs = steps_obs.run_Calc_Blocks(ibls_obs,gibls_obs,lons_obs,daynum_obs)

    if "CALCBLOCKS" in steps_list_fcst:
        if "GIBL" not in steps_list_fcst:
            raise ValueError('Must run Forecast GIBLs before calculating blocks.')
        print('Computing Forecast Blocks')
        block_freq_fcst = steps_fcst.run_Calc_Blocks(ibls_fcst,gibls_fcst,lons_fcst,daynum_fcst)

    # Write out a Blocking MPR file if both obs and forecast blocking calculation performed
    if "CALCBLOCKS" in steps_list_obs and "CALCBLOCKS" in steps_list_fcst:
        b_mpr_outdir = os.path.join(mpr_dir,'Blocks')
        if not os.path.exists(b_mpr_outdir):
            os.makedirs(b_mpr_outdir)
        # Print Blocks to output matched pair file
        modname = os.environ.get('MODEL_NAME','GFS')
        maskname = os.environ.get('MASK_NAME','FULL')
        blocks_outfile_prefix = os.path.join(b_mpr_outdir,'blocking_stat_'+modname)
        cbls_avg = np.nanmean(cbls_obs,axis=0)
        write_mpr_file(block_freq_obs,block_freq_fcst,cbls_avg,lons_obs,ibl_time_obs,ibl_time_fcst,modname,
            'NA','Blocks','block','Z500','Blocks','block','Z500',maskname,'500',blocks_outfile_prefix)


    # Plot Blocking Frequency
    if "PLOTBLOCKS" in steps_list_obs:
        if "CALCBLOCKS" not in steps_list_obs:
            raise ValueError('Must compute observed blocks before plotting them.')
        print('Plotting Obs Blocks')
        blocking_plot_title = os.environ.get('OBS_BLOCKING_PLOT_TITLE','Obs Blocking Frequency')
        blocking_plot_outname = os.path.join(oplot_dir,os.environ.get('OBS_BLOCKING_PLOT_OUTPUT_NAME','obs_Block_Freq'))
        pb.plot_blocks(block_freq_obs,gibls_obs,ibls_obs,lons_obs,blocking_plot_title,blocking_plot_outname)
    if "PLOTBLOCKS" in steps_list_fcst:
        if "CALCBLOCKS" not in steps_list_fcst:
            raise ValueError('Must compute forecast blocks before plotting them.')
        print('Plotting Forecast Blocks')
        blocking_plot_title = os.environ.get('FCST_BLOCKING_PLOT_TITLE','Forecast Blocking Frequency')
        blocking_plot_outname = os.path.join(oplot_dir,os.environ.get('FCST_BLOCKING_PLOT_OUTPUT_NAME','fcst_Block_Freq'))
        pb.plot_blocks(block_freq_fcst,gibls_fcst,ibls_fcst,lons_fcst,blocking_plot_title,blocking_plot_outname)


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_mid_lat/UserScript_obsERA_obsOnly_Blocking.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_mid_lat/Blocking. There should be 3 different graphics output to the plot directory in the above location. Each will have png and pdf versions to make for 6 output plots:

* ERA_CBL_avg.png
* ERA_CBL_avg.pdf
* ERA_IBL_Freq_DJF.png
* ERA_IBL_Freq_DJF.pdf
* obs_Block_Freq_DJF.png
* obs_Block_Freq_DJF.pdf

If the pre-processing steps are turned on, regridded data, daily averaged files, running mean files, and anomaly files will also be output to Regrid, Daily,Rmean5d, and Anomaly directories in the ERA directory.

Keywords

Note

  • RegridDataPlaneToolUseCase

  • PCPCombineToolUseCase

  • StatAnalysisToolUseCase

  • S2SAppUseCase

  • S2SMidLatAppUseCase

  • UserScriptUseCase

  • NetCDFFileUseCase

  • GRIB2FileUseCase

  • METcalcpyUseCase

  • METplotpyUseCase

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

Gallery generated by Sphinx-Gallery