Examples

The following examples illustrate practical uses of eof2. Many of the examples have a version for each of

Geopotential height (NAO) [cdms2 interface]

"""
Compute and plot the leading EOF of geopotential height on the 500 hPa
pressure surface over the European/Atlantic sector during winter time.

This example uses the metadata-retaining cdms2 interface.

"""
import cdms2
import cdutil
import matplotlib.pyplot as plt
from mpl_toolkits.basemap import Basemap
import numpy as np

from eof2 import Eof


# Read geopotential height data using the cdms2 module from CDAT. The file
# contains December-February averages of geopotential height at 500 hPa for
# the European/Atlantic domain (80W-40E, 20-90N).
ncin = cdms2.open('../../example_data/hgt_djf.nc', 'r')
z_djf = ncin('z')
ncin.close()

# Compute anomalies by removing the time-mean.
z_djf_mean = cdutil.averager(z_djf, axis='t')
z_djf = z_djf - z_djf_mean
z_djf.id = 'z'

# Create an EOF solver to do the EOF analysis. Square-root of cosine of
# latitude weights are applied before the computation of EOFs.
solver = Eof(z_djf, weights='coslat')

# Retrieve the leading EOF, expressed as the covariance between the leading PC
# time series and the input SLP anomalies at each grid point.
eof1 = solver.eofsAsCovariance(neofs=1)

# Plot the leading EOF expressed as covariance in the European/Atlantic domain.
m = Basemap(projection='ortho', lat_0=60., lon_0=-20.)
lons, lats = eof1.getLongitude()[:], eof1.getLatitude()[:]
x, y = m(*np.meshgrid(lons, lats))
m.contourf(x, y, eof1(squeeze=True), cmap=plt.cm.RdBu_r)
m.drawcoastlines()
m.drawparallels(np.arange(-80, 90, 20))
m.drawmeridians(np.arange(0, 360, 20))
plt.title('EOF1 expressed as covariance', fontsize=16)

plt.show()

(Source code, png, hires.png, pdf)

_images/hgt_example.png

Sea surface temperature (El Niño) [cdms2 interface]

"""
Compute and plot the leading EOF of sea surface temperature in the
central and northern Pacific during winter time.

The spatial pattern of this EOF is the canonical El Nino pattern, and
the associated time series shows large peaks and troughs for well-known
El Nino and La Nina events.

This example uses the metadata-retaining cdms2 interface.

"""
import cdms2
import matplotlib.pyplot as plt
from mpl_toolkits.basemap import Basemap
import numpy as np

from eof2 import Eof


# Read SST anomalies using the cdms2 module from CDAT. The file contains
# November-March averages of SST anomaly in the central and northern Pacific.
ncin = cdms2.open('../../example_data/sst_ndjfm_anom.nc')
sst = ncin('sst')
ncin.close()

# Create an EOF solver to do the EOF analysis. Square-root of cosine of
# latitude weights are applied before the computation of EOFs.
solver = Eof(sst, weights='coslat')

# Retrieve the leading EOF, expressed as the correlation between the leading
# PC time series and the input SST anomalies at each grid point, and the
# leading PC time series itself.
eof1 = solver.eofsAsCorrelation(neofs=1)
pc1 = solver.pcs(npcs=1, pcscaling=1)

# Plot the leading EOF expressed as correlation in the Pacific domain.
m = Basemap(projection='cyl', llcrnrlon=120, llcrnrlat=-20,
        urcrnrlon=260, urcrnrlat=60)
lons, lats = eof1.getLongitude()[:], eof1.getLatitude()[:]
x, y = m(*np.meshgrid(lons, lats))
clevs = np.linspace(-1, 1, 11)
m.contourf(x, y, eof1(squeeze=True), clevs, cmap=plt.cm.RdBu_r)
m.drawcoastlines()
m.drawparallels([-20, 0, 20, 40, 60])
m.drawmeridians([120, 140, 160, 180, 200, 220, 240, 260])
cb = plt.colorbar(orientation='horizontal')
cb.set_label('correlation coefficient', fontsize=12)
plt.title('EOF1 expressed as correlation', fontsize=16)

# Plot the leading PC time series.
plt.figure()
years = range(1962, 2012)
plt.plot(years, pc1, color='b', linewidth=2)
plt.axhline(0, color='k')
plt.title('PC1 Time Series')
plt.xlabel('Year')
plt.ylabel('Normalized Units')
plt.xlim(1962, 2012)
plt.ylim(-3, 3)

plt.show()

(Source code)

Geopotential height (NAO) [numpy interface]

"""
Compute and plot the leading EOF of geopotential height on the 500 hPa
pressure surface over the European/Atlantic sector during winter time.

This example uses the plain numpy interface.

"""
from netCDF4 import Dataset
import matplotlib.pyplot as plt
from mpl_toolkits.basemap import Basemap
import numpy as np

from eof2 import EofSolver


# Read geopotential height data using the netCDF4 module. The file contains
# December-February averages of geopotential height at 500 hPa for the
# European/Atlantic domain (80W-40E, 20-90N).
ncin = Dataset('../../example_data/hgt_djf.nc', 'r')
z_djf = ncin.variables['z'][:]
lons = ncin.variables['longitude'][:]
lats = ncin.variables['latitude'][:]
ncin.close()

# Compute anomalies by removing the time-mean.
z_djf_mean = z_djf.mean(axis=0)
z_djf = z_djf - z_djf_mean

# Create an EOF solver to do the EOF analysis. Square-root of cosine of
# latitude weights are applied before the computation of EOFs.
coslat = np.cos(np.deg2rad(lats)).clip(0.,1.)
wgts = np.sqrt(coslat)[..., np.newaxis]
solver = EofSolver(z_djf, weights=wgts)

# Retrieve the leading EOF, expressed as the covariance between the leading PC
# time series and the input SLP anomalies at each grid point.
eof1 = solver.eofsAsCovariance(neofs=1)

# Plot the leading EOF expressed as covariance in the European/Atlantic domain.
m = Basemap(projection='ortho', lat_0=60., lon_0=-20.)
x, y = m(*np.meshgrid(lons, lats))
m.contourf(x, y, eof1.squeeze(), cmap=plt.cm.RdBu_r)
m.drawcoastlines()
m.drawparallels(np.arange(-80, 90, 20))
m.drawmeridians(np.arange(0, 360, 20))
plt.title('EOF1 expressed as covariance', fontsize=16)

plt.show()

(Source code, png, hires.png, pdf)

_images/hgt_example1.png

Sea surface temperature (El Niño) [numpy interface]

"""
Compute and plot the leading EOF of sea surface temperature in the
central and northern Pacific during winter time.

The spatial pattern of this EOF is the canonical El Nino pattern, and
the associated time series shows large peaks and troughs for well-known
El Nino and La Nina events.

This example uses the plain numpy interface.

"""
from netCDF4 import Dataset
import matplotlib.pyplot as plt
from mpl_toolkits.basemap import Basemap
import numpy as np

from eof2 import EofSolver


# Read SST anomalies using the netCDF4 module. The file contains
# November-March averages of SST anomaly in the central and northern Pacific.
ncin = Dataset('../../example_data/sst_ndjfm_anom.nc', 'r')
sst = ncin.variables['sst'][:]
lons = ncin.variables['longitude'][:]
lats = ncin.variables['latitude'][:]
ncin.close()

# Create an EOF solver to do the EOF analysis. Square-root of cosine of
# latitude weights are applied before the computation of EOFs.
coslat = np.cos(np.deg2rad(lats))
wgts = np.sqrt(coslat)[..., np.newaxis]
solver = EofSolver(sst, weights=wgts)

# Retrieve the leading EOF, expressed as the correlation between the leading
# PC time series and the input SST anomalies at each grid point, and the
# leading PC time series itself.
eof1 = solver.eofsAsCorrelation(neofs=1)
pc1 = solver.pcs(npcs=1, pcscaling=1)

# Plot the leading EOF expressed as correlation in the Pacific domain.
m = Basemap(projection='cyl', llcrnrlon=120, llcrnrlat=-20,
        urcrnrlon=260, urcrnrlat=60)
x, y = m(*np.meshgrid(lons, lats))
clevs = np.linspace(-1, 1, 11)
m.contourf(x, y, eof1.squeeze(), clevs, cmap=plt.cm.RdBu_r)
m.drawcoastlines()
m.drawparallels([-20, 0, 20, 40, 60])
m.drawmeridians([120, 140, 160, 180, 200, 220, 240, 260])
cb = plt.colorbar(orientation='horizontal')
cb.set_label('correlation coefficient', fontsize=12)
plt.title('EOF1 expressed as correlation', fontsize=16)

# Plot the leading PC time series.
plt.figure()
years = range(1962, 2012)
plt.plot(years, pc1, color='b', linewidth=2)
plt.axhline(0, color='k')
plt.title('PC1 Time Series')
plt.xlabel('Year')
plt.ylabel('Normalized Units')
plt.xlim(1962, 2012)
plt.ylim(-3, 3)

plt.show()

(Source code)