xarray.DataArray.xvec.set_crs#

DataArray.xvec.set_crs(variable_crs=None, allow_override=False, **variable_crs_kwargs)#

Set the Coordinate Reference System (CRS) of coordinates backed by GeometryIndex.

Parameters:
variable_crsdict-like or None, optional

A dict where the keys are the names of the coordinates and values target CRS in any format accepted by pyproj.CRS.from_user_input() such as an authority string (e.g. "EPSG:4326"), EPSG code (e.g. 4326) or a WKT string.

allow_overridebool, default False

If the GeometryIndex already has a CRS, allow to replace the existing CRS, even when both are not equal.

**variable_crs_kwargsoptional

The keyword arguments form of variable_crs. One of variable_crs or variable_crs_kwargs must be provided.

Returns:
assignedsame type as caller

A new object with the assigned target CRS.

See also

to_crs

Notes

The underlying geometries are not transformed to this CRS. To transform the geometries to a new CRS, use the to_crs() method.

Examples

The method is the most useful, when there is no CRS assigned (illustrated on a xarray.Dataset but the same is applicable on a xarray.DataArray).

>>> ds = (
...     xr.Dataset(coords={"geom": [shapely.Point(1, 2), shapely.Point(3, 4)]})
...     .xvec.set_geom_indexes("geom")
... )
>>> ds
<xarray.Dataset>
Dimensions:  (geom: 2)
Coordinates:
  * geom     (geom) object POINT (1 2) POINT (3 4)
Data variables:
    *empty*
Indexes:
    geom     GeometryIndex (crs=None)
>>> ds.xvec.set_crs(geom=4326)
<xarray.Dataset>
Dimensions:  (geom: 2)
Coordinates:
  * geom     (geom) object POINT (1 2) POINT (3 4)
Data variables:
    *empty*
Indexes:
    geom     GeometryIndex (crs=EPSG:4326)

It can also be used to overwrite the existing CRS. Note, that in most cases you probably want to use the to_crs() instead is such a case.

>>> ds = (
...     xr.Dataset(coords={"geom": [shapely.Point(1, 2), shapely.Point(3, 4)]})
...     .xvec.set_geom_indexes("geom", crs=4326)
... )
>>> ds
<xarray.Dataset>
Dimensions:  (geom: 2)
Coordinates:
* geom     (geom) object POINT (1 2) POINT (3 4)
Data variables:
    *empty*
Indexes:
    geom     GeometryIndex (crs=EPSG:4326)
>>> ds.xvec.set_crs(geom=3857, allow_override=True)
<xarray.Dataset>
Dimensions:  (geom: 2)
Coordinates:
* geom     (geom) object POINT (1 2) POINT (3 4)
Data variables:
    *empty*
Indexes:
    geom     GeometryIndex (crs=EPSG:3857)

See that only the CRS has changed, not the geometries.