Read and Write netCDF Files#

Reading#

Read data from a NetCDF file into a RasterModelGrid.

Read netcdf#

read_netcdf(nc_file[, grid, name, ...])

Create a RasterModelGrid from a netcdf file.

read_netcdf(nc_file, grid=None, name=None, just_grid=False, halo=0, nodata_value=-9999.0)[source]#

Create a RasterModelGrid from a netcdf file.

Create a new RasterModelGrid from the netcdf file, nc_file. If the netcdf file also contains data, add that data to the grid’s fields. To create a new grid without any associated data from the netcdf file, set the just_grid keyword to True.

A halo can be added with the keyword halo.

If you want the fields to be added to an existing grid, it can be passed to the keyword argument grid.

Parameters:
  • nc_file (str) – Name of a netcdf file.

  • grid (grid , optional) – Adds data to an existing grid instead of creating a new one.

  • name (str, optional) – Add only fields with NetCDF variable name to the grid. Default is to add all NetCDF varibles to the grid.

  • just_grid (boolean, optional) – Create a new grid but don’t add value data.

  • halo (integer, optional) – Adds outer border of depth halo to the grid.

  • nodata_value (float, optional) – Value that indicates an invalid value. Default is -9999.

Returns:

A newly-created RasterModelGrid.

Return type:

RasterModelGrid

Examples

Import read_netcdf and the path to an example netcdf file included with landlab.

>>> from landlab.io.netcdf import read_netcdf

Create a new grid from the netcdf file. The example grid is a uniform rectilinear grid with 4 rows and 3 columns of nodes with unit spacing. The data file also contains data defined at the nodes for the grid for a variable called, surface__elevation.

>>> grid = read_netcdf("test-netcdf4.nc")  
>>> grid.shape == (4, 3)  
True
>>> grid.dy, grid.dx  
(1.0, 1.0)
>>> list(grid.at_node.keys())  
['surface__elevation']
>>> grid.at_node["surface__elevation"]  
array([  0.,   1.,   2.,   3.,   4.,   5.,   6.,   7.,   8.,   9.,  10.,
        11.])

read_netcdf will try to determine the format of the netcdf file. For example, the same call will also work for netcdf3-formatted files.

>>> grid = read_netcdf("test-netcdf3-64bit.nc")  
>>> grid.shape == (4, 3)  
True
>>> grid.dy, grid.dx  
(1.0, 1.0)

A more complicated example might add data with a halo to an existing grid. Note that the lower left corner must be specified correctly for the data and the grid to align correctly.

>>> from landlab import RasterModelGrid
>>> grid = RasterModelGrid((6, 5), xy_of_lower_left=(-1.0, -1.0))  
>>> grid = read_netcdf(
...     "test-netcdf4.nc",
...     grid=grid,
...     halo=1,
...     nodata_value=-1,
... )  
>>> grid.at_node["surface__elevation"].reshape(grid.shape)  
array([[ -1.,  -1.,  -1.,  -1.,  -1.],
       [ -1.,   0.,   1.,   2.,  -1.],
       [ -1.,   3.,   4.,   5.,  -1.],
       [ -1.,   6.,   7.,   8.,  -1.],
       [ -1.,   9.,  10.,  11.,  -1.],
       [ -1.,  -1.,  -1.,  -1.,  -1.]])

Writing#

Write structured grids to NetCDF files.

Write netcdf#

write_netcdf(path, grid[, attrs, append, ...])

Write landlab fields to netcdf.

write_netcdf(path, grid, attrs=None, append=False, format='NETCDF3_64BIT', names=None, at=None, time=None, raster=False)[source]#

Write landlab fields to netcdf.

Write the data and grid information for grid to path as NetCDF. If the append keyword argument in True, append the data to an existing file, if it exists. Otherwise, clobber an existing files.

Parameters:
  • path (str) – Path to output file.

  • grid (RasterModelGrid) – Landlab RasterModelGrid object that holds a grid and associated values.

  • append (boolean, optional) – Append data to an existing file, otherwise clobber the file.

  • format ({'NETCDF3_CLASSIC', 'NETCDF3_64BIT', 'NETCDF4_CLASSIC', 'NETCDF4'}) – Format of output netcdf file.

  • attrs (dict) – Attributes to add to netcdf file.

  • names (iterable of str, optional) – Names of the fields to include in the netcdf file. If not provided, write all fields.

  • at ({'node', 'cell'}, optional) – The location where values are defined.

  • raster (bool, optional) – Indicate whether spatial dimensions are written as full value arrays (default) or just as coordinate dimensions.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> from landlab.io.netcdf import write_netcdf

Create a uniform rectilinear grid with four rows and 3 columns, and add some data fields to it.

>>> rmg = RasterModelGrid((4, 3))
>>> rmg.at_node["topographic__elevation"] = np.arange(12.0)
>>> rmg.at_node["uplift_rate"] = 2.0 * np.arange(12.0)

Create a temporary directory to write the netcdf file into.

>>> import tempfile, os
>>> temp_dir = tempfile.mkdtemp()
>>> os.chdir(temp_dir)

Write the grid to a netcdf3 file but only include the uplift_rate data in the file.

>>> write_netcdf("test.nc", rmg, format="NETCDF3_64BIT", names="uplift_rate")

Read the file back in and check its contents.

>>> from scipy.io import netcdf
>>> fp = netcdf.netcdf_file("test.nc", "r")
>>> "uplift_rate" in fp.variables
True
>>> "topographic__elevation" in fp.variables
False
>>> fp.variables["uplift_rate"][:].flatten().astype("=f8")
array([  0.,   2.,   4.,   6.,   8.,  10.,  12.,  14.,  16.,  18.,  20.,
        22.])
>>> rmg.at_cell["air__temperature"] = np.arange(2.0)
>>> write_netcdf(
...     "test-cell.nc",
...     rmg,
...     format="NETCDF3_64BIT",
...     names="air__temperature",
...     at="cell",
... )
write_raster_netcdf(path, fields, attrs=None, append=False, time=None, format='NETCDF4', names=None, at=None)[source]#

Write Raster Model Grid landlab fields to netcdf.

Write the data and grid information for fields to path as NetCDF.

This method is for Raster Grids only and takes advantage of regular x and y spacing to save memory.

Rather that writing x and y of node locations at all (nr x nc) locations, it writes a 1D array each for x and y.

A more modern version of this might write x and y location as a netcdf coordinate. However, the original version of this function wrote x and y as data variables rather than coordinates.

If the append keyword argument in True, append the data to an existing file, if it exists. Otherwise, clobber an existing files.

Parameters:
  • path (str) – Path to output file.

  • fields (field-like) – Landlab field object that holds a grid and associated values. This must be a Raster type.

  • append (boolean, optional) – Append data to an existing file, otherwise clobber the file.

  • time (float, optional) – Add a time to the time variable.

  • format ({'NETCDF4'}) – Format of output netcdf file.

  • attrs (dict) – Attributes to add to netcdf file.

  • names (iterable of str, optional) – Names of the fields to include in the netcdf file. If not provided, write all fields.

  • at ({'node'}, optional) – The location where values are defined. Presently only implemented for type ‘node’.

Examples

>>> import numpy as np
>>> from landlab import RasterModelGrid
>>> from landlab.io.netcdf import write_raster_netcdf

Create a uniform rectilinear grid with four rows and 3 columns, and add some data fields to it.

>>> rmg = RasterModelGrid((4, 3))
>>> rmg.shape
(4, 3)
>>> rmg.at_node["topographic__elevation"] = np.arange(12.0)
>>> rmg.at_node["uplift_rate"] = 2.0 * np.arange(12.0)

Create a temporary directory to write the netcdf file into.

>>> import tempfile, os
>>> temp_dir = tempfile.mkdtemp()
>>> os.chdir(temp_dir)

Write the grid to a netcdf4 file but only include the uplift_rate data in the file.

>>> write_raster_netcdf(
...     "test.nc",
...     rmg,
...     format="NETCDF3_64BIT",
...     names="uplift_rate",
... )

Read the file back in and check its contents.

>>> from scipy.io import netcdf
>>> fp = netcdf.netcdf_file("test.nc", "r")
>>> "uplift_rate" in fp.variables
True
>>> "topographic__elevation" in fp.variables
False
>>> fp.variables["uplift_rate"][:].flatten().astype("=f8")
array([  0.,   2.,   4.,   6.,   8.,  10.,  12.,  14.,  16.,  18.,  20.,
        22.])
>>> fp.variables["x"][:].astype("=f8")
array([ 0.,  1.,  2.])
>>> fp.variables["y"][:].astype("=f8")
array([ 0.,  1.,  2.,  3.])

Read now with read_netcdf

>>> from landlab.io.netcdf import read_netcdf
>>> grid = read_netcdf("test.nc")
>>> grid.shape
(4, 3)
>>> grid.x_of_node
array([ 0.,  1.,  2.,  0.,  1.,  2.,  0.,  1.,  2.,  0.,  1.,  2.])
>>> grid.y_of_node
array([ 0.,  0.,  0.,  1.,  1.,  1.,  2.,  2.,  2.,  3.,  3.,  3.])
>>> grid.at_node["uplift_rate"]
array([  0.,   2.,   4.,   6.,   8.,  10.,  12.,  14.,  16.,  18.,  20.,
        22.])

Error Handling#

Exceptions to raise for the netcdf module.

exception Error[source]#

Bases: Exception

Base class for errors in this package.

exception NotRasterGridError[source]#

Bases: Error

Raise if grid is not uniform rectilinear.

Raise this error if the grid defined in the netcdf file is not uniform rectilinear with constant spacing in all dimensions.