Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
20 changes: 20 additions & 0 deletions docs/make_projection.py
Original file line number Diff line number Diff line change
Expand Up @@ -41,6 +41,25 @@ def plate_carree_plot():
ax.gridlines()


def igh_plot():
import matplotlib.pyplot as plt
import cartopy.crs as ccrs

fig = plt.figure(figsize=(6.9228, 6))

ax1 = fig.add_subplot(2, 1, 1,
projection=ccrs.InterruptedGoodeHomolosine(
emphasis='land'))
ax1.coastlines(resolution='110m')
ax1.gridlines()

ax2 = fig.add_subplot(2, 1, 2,
projection=ccrs.InterruptedGoodeHomolosine(
central_longitude=-160, emphasis='ocean'))
ax2.coastlines(resolution='110m')
ax2.gridlines()


def utm_plot():
import matplotlib.pyplot as plt
import cartopy.crs as ccrs
Expand All @@ -59,6 +78,7 @@ def utm_plot():

MULTI_PLOT_CASES = {
ccrs.PlateCarree: plate_carree_plot,
ccrs.InterruptedGoodeHomolosine: igh_plot,
ccrs.UTM: utm_plot,
}

Expand Down
2 changes: 1 addition & 1 deletion docs/source/contributors.rst
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ the package wouldn't be as rich or diverse as it is today:
* Sadie Bartholomew
* Kacper Makuch
* Stephane Raynaud

* John Krasting

Thank you!

Expand Down
17 changes: 13 additions & 4 deletions docs/source/reference/projections.rst
Original file line number Diff line number Diff line change
Expand Up @@ -271,10 +271,19 @@ InterruptedGoodeHomolosine
import matplotlib.pyplot as plt
import cartopy.crs as ccrs

plt.figure(figsize=(6.9228, 3))
ax = plt.axes(projection=ccrs.InterruptedGoodeHomolosine())
ax.coastlines(resolution='110m')
ax.gridlines()
fig = plt.figure(figsize=(6.9228, 6))

ax1 = fig.add_subplot(2, 1, 1,
projection=ccrs.InterruptedGoodeHomolosine(
emphasis='land'))
ax1.coastlines(resolution='110m')
ax1.gridlines()

ax2 = fig.add_subplot(2, 1, 2,
projection=ccrs.InterruptedGoodeHomolosine(
central_longitude=-160, emphasis='ocean'))
ax2.coastlines(resolution='110m')
ax2.gridlines()


RotatedPole
Expand Down
52 changes: 47 additions & 5 deletions lib/cartopy/crs.py
Original file line number Diff line number Diff line change
Expand Up @@ -1930,17 +1930,59 @@ def transform_points(self, src_crs, x, y, z=None):


class InterruptedGoodeHomolosine(Projection):
def __init__(self, central_longitude=0, globe=None):
proj4_params = [('proj', 'igh'), ('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)
"""
Composite equal-area projection empahsizing either land or
ocean features.

Original Reference:
Goode, J. P., 1925: The Homolosine Projection: A new device for
portraying the Earth's surface entire. Annals of the
Association of American Geographers, 15:3, 119-125,
DOI: 10.1080/00045602509356949

A central_longitude value of -160 is recommended for the oceanic view.

"""
def __init__(self, central_longitude=0, globe=None, emphasis='land'):
"""
Parameters
----------
central_longitude : float, optional
The central longitude, by default 0
globe : :class:`cartopy.crs.Globe`, optional
If omitted, a default Globe object is created, by default None
emphasis : str, optional
Options 'land' and 'ocean' are available, by default 'land'
"""

if emphasis == 'land':
proj4_params = [('proj', 'igh'), ('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)

elif emphasis == 'ocean':
if PROJ4_VERSION < (7, 1, 0):
_proj_ver = '.'.join(str(v) for v in PROJ4_VERSION)
raise ValueError('The Interrupted Goode Homolosine ocean '
'projection requires Proj version 7.1.0, '
'but you are using ' + _proj_ver)
proj4_params = [('proj', 'igh_o'), ('lon_0', central_longitude)]
super().__init__(proj4_params, globe=globe)

else:
msg = '`emphasis` needs to be either \'land\' or \'ocean\''
raise ValueError(msg)

minlon, maxlon = self._determine_longitude_bounds(central_longitude)
epsilon = 1e-10

# Obtain boundary points
n = 31
top_interrupted_lons = (-40.0,)
bottom_interrupted_lons = (80.0, -20.0, -100.0)
if emphasis == 'land':
top_interrupted_lons = (-40.0,)
bottom_interrupted_lons = (80.0, -20.0, -100.0)
elif emphasis == 'ocean':
top_interrupted_lons = (-90.0, 60.0)
bottom_interrupted_lons = (90.0, -60.0)
lons = np.empty(
(2 + 2 * len(top_interrupted_lons + bottom_interrupted_lons)) * n +
1)
Expand Down
93 changes: 58 additions & 35 deletions lib/cartopy/tests/crs/test_interrupted_goode_homolosine.py
Original file line number Diff line number Diff line change
Expand Up @@ -16,38 +16,61 @@
from .helpers import check_proj_params


def test_default():
igh = ccrs.InterruptedGoodeHomolosine()
other_args = {'ellps=WGS84', 'lon_0=0'}
check_proj_params('igh', igh, other_args)

assert_almost_equal(np.array(igh.x_limits),
[-20037508.3427892, 20037508.3427892])
assert_almost_equal(np.array(igh.y_limits),
[-8683259.7164347, 8683259.7164347])


def test_eccentric_globe():
globe = ccrs.Globe(semimajor_axis=1000, semiminor_axis=500,
ellipse=None)
igh = ccrs.InterruptedGoodeHomolosine(globe=globe)
other_args = {'a=1000', 'b=500', 'lon_0=0'}
check_proj_params('igh', igh, other_args)

assert_almost_equal(np.array(igh.x_limits),
[-3141.5926536, 3141.5926536])
assert_almost_equal(np.array(igh.y_limits),
[-1361.410035, 1361.410035])


@pytest.mark.parametrize('lon', [-10.0, 10.0])
def test_central_longitude(lon):
igh = ccrs.InterruptedGoodeHomolosine(central_longitude=lon)
other_args = {'ellps=WGS84', 'lon_0={}'.format(lon)}
check_proj_params('igh', igh, other_args)

assert_almost_equal(np.array(igh.x_limits),
[-20037508.3427892, 20037508.3427892],
decimal=5)
assert_almost_equal(np.array(igh.y_limits),
[-8683259.7164347, 8683259.7164347])
@pytest.mark.parametrize("emphasis", ["land", "ocean"])
def test_default(emphasis):
if emphasis == "ocean" and ccrs.PROJ4_VERSION < (7, 1, 0):
pytest.skip()
igh = ccrs.InterruptedGoodeHomolosine(emphasis=emphasis)
other_args = {"ellps=WGS84", "lon_0=0"}
if emphasis == "land":
check_proj_params("igh", igh, other_args)
elif emphasis == "ocean":
check_proj_params("igh_o", igh, other_args)
assert_almost_equal(
np.array(igh.x_limits), [-20037508.3427892, 20037508.3427892]
)
assert_almost_equal(
np.array(igh.y_limits), [-8683259.7164347, 8683259.7164347]
)


@pytest.mark.parametrize("emphasis", ["land", "ocean"])
def test_eccentric_globe(emphasis):
if emphasis == "ocean" and ccrs.PROJ4_VERSION < (7, 1, 0):
pytest.skip()
globe = ccrs.Globe(semimajor_axis=1000, semiminor_axis=500, ellipse=None)
igh = ccrs.InterruptedGoodeHomolosine(globe=globe, emphasis=emphasis)
other_args = {"a=1000", "b=500", "lon_0=0"}
if emphasis == "land":
check_proj_params("igh", igh, other_args)
elif emphasis == "ocean":
check_proj_params("igh_o", igh, other_args)

assert_almost_equal(np.array(igh.x_limits), [-3141.5926536, 3141.5926536])
assert_almost_equal(np.array(igh.y_limits), [-1361.410035, 1361.410035])


@pytest.mark.parametrize(
("emphasis", "lon"),
[("land", -10.0), ("land", 10.0), ("ocean", -10.0), ("ocean", 10.0)],
)
def test_central_longitude(emphasis, lon):
if emphasis == "ocean" and ccrs.PROJ4_VERSION < (7, 1, 0):
pytest.skip()
igh = ccrs.InterruptedGoodeHomolosine(
central_longitude=lon, emphasis=emphasis
)
other_args = {"ellps=WGS84", "lon_0={}".format(lon)}
if emphasis == "land":
check_proj_params("igh", igh, other_args)
elif emphasis == "ocean":
check_proj_params("igh_o", igh, other_args)

assert_almost_equal(
np.array(igh.x_limits),
[-20037508.3427892, 20037508.3427892],
decimal=5,
)
assert_almost_equal(
np.array(igh.y_limits), [-8683259.7164347, 8683259.7164347]
)
Loading
Sorry, something went wrong. Reload?
Sorry, we cannot display this file.
Sorry, this file is invalid so it cannot be displayed.
Loading
Sorry, something went wrong. Reload?
Sorry, we cannot display this file.
Sorry, this file is invalid so it cannot be displayed.
21 changes: 21 additions & 0 deletions lib/cartopy/tests/mpl/test_crs.py
Original file line number Diff line number Diff line change
Expand Up @@ -12,6 +12,27 @@
from cartopy.tests.mpl import ImageTesting


@pytest.mark.natural_earth
@ImageTesting(["igh_land"])
def test_igh_land():
crs = ccrs.InterruptedGoodeHomolosine(emphasis="land")
ax = plt.axes(projection=crs)
ax.coastlines()
ax.gridlines()


@pytest.mark.skipif(ccrs.PROJ4_VERSION < (7, 1, 0), reason="Proj is too old.")
@pytest.mark.natural_earth
@ImageTesting(["igh_ocean"])
def test_igh_ocean():
crs = ccrs.InterruptedGoodeHomolosine(
central_longitude=-160, emphasis="ocean"
)
ax = plt.axes(projection=crs)
ax.coastlines()
ax.gridlines()


@pytest.mark.natural_earth
@ImageTesting(['lambert_conformal_south'])
def test_lambert_south():
Expand Down
2 changes: 2 additions & 0 deletions lib/cartopy/tests/mpl/test_examples.py
Original file line number Diff line number Diff line change
Expand Up @@ -48,6 +48,8 @@ def test_global_map():


contour_image = 'contour_label' if MPL_VERSION < '3.4' else 'contour_label_3.4'


@pytest.mark.natural_earth
@ExampleImageTesting([contour_image], tolerance=0)
def test_contour_label():
Expand Down
2 changes: 1 addition & 1 deletion lib/cartopy/tests/mpl/test_gridliner.py
Original file line number Diff line number Diff line change
Expand Up @@ -32,7 +32,7 @@
ccrs.Robinson,
ccrs.Sinusoidal,
ccrs.Stereographic,
ccrs.InterruptedGoodeHomolosine,
(ccrs.InterruptedGoodeHomolosine, dict(emphasis='land')),
(ccrs.RotatedPole,
dict(pole_longitude=180.0,
pole_latitude=36.0,
Expand Down
2 changes: 1 addition & 1 deletion lib/cartopy/tests/mpl/test_img_transform.py
Original file line number Diff line number Diff line change
Expand Up @@ -96,7 +96,7 @@ def test_regrid_image():
# Target grid
target_nx = 300
target_ny = 300
target_proj = ccrs.InterruptedGoodeHomolosine()
target_proj = ccrs.InterruptedGoodeHomolosine(emphasis='land')
target_x, target_y, target_extent = im_trans.mesh_projection(target_proj,
target_nx,
target_ny)
Expand Down
4 changes: 2 additions & 2 deletions lib/cartopy/tests/mpl/test_mpl_integration.py
Original file line number Diff line number Diff line change
Expand Up @@ -212,7 +212,7 @@ def test_multiple_projections():
ccrs.SouthPolarStereo(),
ccrs.Orthographic(),
ccrs.Mollweide(),
ccrs.InterruptedGoodeHomolosine(),
ccrs.InterruptedGoodeHomolosine(emphasis='land'),
ccrs.EckertI(),
ccrs.EckertII(),
ccrs.EckertIII(),
Expand Down Expand Up @@ -644,7 +644,7 @@ def test_pcolormesh_goode_wrap():
X, Y = np.meshgrid(*[np.deg2rad(c) for c in (x, y)])
Z = np.cos(Y) + 0.375 * np.sin(2. * X)
Z = Z[:-1, :-1]
ax = plt.axes(projection=ccrs.InterruptedGoodeHomolosine())
ax = plt.axes(projection=ccrs.InterruptedGoodeHomolosine(emphasis='land'))
ax.coastlines()
ax.pcolormesh(x, y, Z, transform=ccrs.PlateCarree())

Expand Down