Skip to content

Commit

Permalink
Support different BSF offsets for different plots
Browse files Browse the repository at this point in the history
This requires a new section to be defined for each plot that doesn't
use the default offset.  This new section also allows (in fact
requires) a new colorbar/colormap to be defined.
  • Loading branch information
xylar committed Mar 5, 2025
1 parent c3c60a4 commit bc10aae
Show file tree
Hide file tree
Showing 3 changed files with 149 additions and 30 deletions.
12 changes: 7 additions & 5 deletions mpas_analysis/default.cfg
Original file line number Diff line number Diff line change
Expand Up @@ -265,7 +265,7 @@ comparisonSubpolarNorthAtlanticResolution = 20.
# The comparison FRIS Antarctic polar stereographic grid size and resolution in km
comparisonFrisBounds = [-1800., -400., 100., 1500.]
# comparisonFrisWidth = 6000.
comparisonFrisResolution = 10.
comparisonFrisResolution = 4.

# interpolation order for model and observation results. Likely values are
# 'bilinear', 'neareststod' (nearest neighbor) or 'conserve'
Expand Down Expand Up @@ -1596,12 +1596,14 @@ seasons = ['ANN']
# comparison grid(s) on which to plot analysis
comparisonGrids = ['latlon', 'subpolar_north_atlantic']

# list of tuples(pairs) of depths (min, max) to integrate horizontal transport over
# list of tuples(pairs) of depths (min, max) to integrate horizontal transport
# over
depthRanges = [(10.0, -10000.0), (10.0, -2000.0)]

# minimum latitude (degrees) above which the mean BSF on boundary vertices
# averages to zero
minLatitudeForZeroBSF = -45.0
# minimum and maximum latitude (degrees) between which the mean BSF on boundary
# vertices averages to zero
latitudeRangeForZeroBSF = (-45.0, 90.0)


[climatologyMapOHCAnomaly]
## options related to plotting horizontally remapped climatologies of
Expand Down
94 changes: 70 additions & 24 deletions mpas_analysis/ocean/climatology_map_bsf.py
Original file line number Diff line number Diff line change
Expand Up @@ -17,6 +17,7 @@
from mpas_analysis.shared.climatology import RemapMpasClimatologySubtask
from mpas_analysis.shared.plot import PlotClimatologyMapSubtask
from mpas_analysis.ocean.utility import compute_zmid
from mpas_analysis.shared.projection import comparison_grid_option_suffixes


class ClimatologyMapBSF(AnalysisTask):
Expand Down Expand Up @@ -72,8 +73,6 @@ def __init__(self, config, mpas_climatology_task, control_config=None):
raise ValueError(f'config section {section_name} does not contain '
f'valid list of comparison grids')

mpas_field_name = field_name

for min_depth, max_depth in depth_ranges:
depth_range_string = f'{min_depth:g}_to_{max_depth:g}m'
if max_depth >= -6000. or min_depth < 0.:
Expand All @@ -99,15 +98,29 @@ def __init__(self, config, mpas_climatology_task, control_config=None):
remap_observations_subtask = None
if control_config is None:
ref_title_label = None
ref_field_name = None
diff_title_label = 'Model - Observations'

else:
control_run_name = control_config.get('runs', 'mainRunName')
ref_title_label = f'Control: {control_run_name}'
ref_field_name = mpas_field_name
diff_title_label = 'Main - Control'
for comparison_grid_name in comparison_grid_names:
grid_suffix = \
comparison_grid_option_suffixes[comparison_grid_name]
config_section_name = f'{self.taskName}{grid_suffix}'
if config.has_section(config_section_name):
# if this comparison grid has its own section, there is a
# version of the BSF that has been offset for this region
# and an associated colorbar/colormap
mpas_field_name = f'{field_name}{grid_suffix}'
else:
config_section_name = self.taskName
mpas_field_name = field_name
if control_config is None:
ref_field_name = None
else:
ref_field_name = mpas_field_name

for season in seasons:
# make a new subtask for this season and comparison grid
subtask_name = f'plot{season}_{comparison_grid_name}' \
Expand All @@ -129,7 +142,8 @@ def __init__(self, config, mpas_climatology_task, control_config=None):
galleryGroup='Horizontal Streamfunction',
groupSubtitle=None,
groupLink='bsf',
galleryName=None)
galleryName=None,
configSectionName=config_section_name)

self.add_subtask(subtask)

Expand Down Expand Up @@ -258,6 +272,7 @@ def customize_masked_climatology(self, climatology, season):
the modified climatology data set
"""
logger = self.logger
config = self.config

ds_mesh = xr.open_dataset(self.restartFileName)
ds_mesh = ds_mesh[['cellsOnEdge', 'cellsOnVertex', 'nEdgesOnCell',
Expand All @@ -276,6 +291,28 @@ def customize_masked_climatology(self, climatology, season):

climatology = climatology.drop_vars(self.variableList)

# offset the BSF for specific comparison grids if defined
for comparison_grid_name in self.comparisonDescriptors.keys():
grid_suffix = \
comparison_grid_option_suffixes[comparison_grid_name]
config_section_name = f'{self.taskName}{grid_suffix}'
if config.has_section(config_section_name):
# if this comparison grid has its own section, there is a
# version of the BSF that has been offset for this region
# and an associated colorbar/colormap
mpas_field_name = \
f'barotropicStreamfunction{grid_suffix}'

lat_range = config.getexpression(
config_section_name, 'latitudeRangeForZeroBSF')
climatology[mpas_field_name] = _shift_bsf(
bsf_vertex, lat_range, ds_mesh.cellsOnVertex - 1,
ds_mesh.latVertex)
climatology[mpas_field_name].attrs['units'] = 'Sv'
climatology[mpas_field_name].attrs['description'] = \
f'barotropic streamfunction at vertices, offset for ' \
f'{grid_suffix} plots'

return climatology

def _compute_vert_integ_velocity(self, ds_mesh, ds):
Expand Down Expand Up @@ -384,8 +421,8 @@ def _compute_barotropic_streamfunction_vertex(self, ds_mesh, ds):
self.logger.info('vertically integrated vorticity computed.')

config = self.config
min_lat = config.getfloat(
'climatologyMapBSF', 'minLatitudeForZeroBSF')
lat_range = config.getexpression(
'climatologyMapBSF', 'latitudeRangeForZeroBSF')

nvertices = ds_mesh.sizes['nVertices']
vertex_degree = ds_mesh.sizes['vertexDegree']
Expand Down Expand Up @@ -468,25 +505,34 @@ def _compute_barotropic_streamfunction_vertex(self, ds_mesh, ds):
bsf_vertex = xr.DataArray(-1e-6 * solution[0:-1],
dims=('nVertices',))

# Now, we actually want the value of the streamfunciton to be
# as close as possible to zero at boundaries north of the ACC
is_boundary_cov = cells_on_vertex == -1
boundary_vertices = is_boundary_cov.sum(dim='vertexDegree') > 0
bsf_vertex = _shift_bsf(bsf_vertex, lat_range, cells_on_vertex,
ds_mesh.latVertex)

boundary_vertices = np.logical_and(
boundary_vertices,
ds_mesh.latVertex > np.deg2rad(min_lat)
)
return bsf_vertex

# convert from boolean mask to indices
boundary_vertices = np.flatnonzero(boundary_vertices.values)

mean_boundary_bsf = bsf_vertex.isel(nVertices=boundary_vertices).mean()
print()
print(f'Mean BSF on boundary north of {min_lat}:')
print(f' {mean_boundary_bsf.values}')
print('Subtracting this mean.')
def _shift_bsf(bsf_vertex, lat_range, cells_on_vertex, lat_vertex):
"""
Shift the barotropic streamfunction to be zero at the boundary over
the given latitude range
"""
is_boundary_cov = cells_on_vertex == -1
boundary_vertices = is_boundary_cov.sum(dim='vertexDegree') > 0

bsf_vertex -= mean_boundary_bsf
boundary_vertices = np.logical_and(
boundary_vertices,
lat_vertex >= np.deg2rad(lat_range[0])
)
boundary_vertices = np.logical_and(
boundary_vertices,
lat_vertex <= np.deg2rad(lat_range[1])
)

return bsf_vertex
# convert from boolean mask to indices
boundary_vertices = np.flatnonzero(boundary_vertices.values)

mean_boundary_bsf = bsf_vertex.isel(nVertices=boundary_vertices).mean()

bsf_shifted = bsf_vertex - mean_boundary_bsf

return bsf_shifted
73 changes: 72 additions & 1 deletion mpas_analysis/polar_regions.cfg
Original file line number Diff line number Diff line change
Expand Up @@ -611,7 +611,78 @@ yLim = [-600., -5.]

# comparison grid(s) on which to plot analysis
comparisonGrids = ['latlon', 'arctic_extended', 'antarctic_extended',
'subpolar_north_atlantic']
'subpolar_north_atlantic', 'fris']

[climatologyMapBSFAntarcticExtended]
## options related to plotting Antarctic climatologies of the barotropic
## streamfunction (BSF)

# colormap for model/observations
colormapNameResult = blue-orange-div
# whether the colormap is indexed or continuous
colormapTypeResult = continuous
# color indices into colormapName for filled contours
# the type of norm used in the colormap
normTypeResult = symLog
# A dictionary with keywords for the norm
normArgsResult = {'linthresh': 30., 'linscale': 0.5, 'vmin': -150., 'vmax': 150.}
colorbarTicksResult = [-100.,-40., -20., -10., 0., 10., 20., 40., 100.]
# Adding contour lines to the figure
contourLevelsResult = np.arange(-150., 150.1, 10.)
contourThicknessResult = 0.5
contourColorResult = black
# Add arrows to contour lines
# whether to include arrows on the contour lines showing the direction of flow
arrowsOnContourResult = True
# colormap for differences
colormapNameDifference = cmo.balance
# whether the colormap is indexed or continuous
colormapTypeDifference = continuous
# the type of norm used in the colormap
normTypeDifference = linear
# A dictionary with keywords for the norm
normArgsDifference = {'vmin': -10., 'vmax': 10.}
# colorbarTicksDifference = numpy.linspace(-10., 10., 9)

# minimum and maximum latitude (degrees) between which the mean BSF on boundary
# vertices averages to zero
latitudeRangeForZeroBSF = (-90.0, -60.0)

[climatologyMapBSFFris]
## options related to plotting Filchner-Ronne climatologies of the barotropic
## streamfunction (BSF)

# colormap for model/observations
colormapNameResult = blue-orange-div
# whether the colormap is indexed or continuous
colormapTypeResult = continuous
# color indices into colormapName for filled contours
# the type of norm used in the colormap
normTypeResult = symLog
# A dictionary with keywords for the norm
normArgsResult = {'linthresh': 1., 'linscale': 0.5, 'vmin': -20., 'vmax': 20.}
colorbarTicksResult = [-20., -10.,-4., -2., -1., 0., 1., 2., 4., 10., 20.]
# Adding contour lines to the figure
contourLevelsResult = np.arange(-5., 5.1, 0.5)
contourThicknessResult = 0.5
contourColorResult = black
# Add arrows to contour lines
# whether to include arrows on the contour lines showing the direction of flow
arrowsOnContourResult = True
# colormap for differences
colormapNameDifference = cmo.balance
# whether the colormap is indexed or continuous
colormapTypeDifference = continuous
# the type of norm used in the colormap
normTypeDifference = linear
# A dictionary with keywords for the norm
normArgsDifference = {'vmin': -2., 'vmax': 2.}
# colorbarTicksDifference = numpy.linspace(-10., 10., 9)

# minimum and maximum latitude (degrees) between which the mean BSF on boundary
# vertices averages to zero
latitudeRangeForZeroBSF = (-90.0, -60.0)


[climatologyMapAntarcticMelt]
## options related to plotting horizontally regridded maps of Antarctic
Expand Down

0 comments on commit bc10aae

Please sign in to comment.