starepandas.STAREDataFrame#
- class STAREDataFrame(*args, sids=None, add_sids=False, level=None, trixels=None, add_trixels=False, n_partitions=1, **kwargs)#
- add_trixels(n_partitions=1, num_workers=None, inplace=False, wrap_lon=True)#
Combination of make_trixels() and set_trixels()
- clear_to_level(inplace=False)#
Clears location bits to level
- Parameters
- inplace: bool
If True, modifies the DataFrame in place (do not create a new object).
Examples
>>> sids = [2299437706637111721, 2299435211084507593, 2299566194809236969] >>> sdf = starepandas.STAREDataFrame(sids=sids) >>> sdf.clear_to_level(inplace=False) sids 0 2299437254470270985 1 2299435055447015433 2 2299564797819093001
- drop_na_sids(inplace=False)#
Drop all rows that have NA values for the SIDs and cast the column to numpy.int64
- hex()#
Returns the hex16 representation of the stare column
Examples
>>> sdf = starepandas.STAREDataFrame(sids=[2251799813685252, 4503599627370500]) >>> sdf.hex() ['0x0008000000000004', '0x0010000000000004']
>>> sdf = starepandas.STAREDataFrame(sids=[[2251799813685252, 4503599627370500], ... [4604930618986332164, 4607182418800017412]]) >>> sdf.hex() [['0x0008000000000004', '0x0010000000000004'], ['0x3fe8000000000004', '0x3ff0000000000004']]
- make_sids(level, convex=False, force_ccw=True, n_partitions=1)#
Generates and returns the STARE representation of each feauture.
- Parameters
- level: int; 0<=level<=27
STARE level to use for the STARE lookup
- convex: bool
Toggle if STARE indices for the convex hull rather than the G-Ring should be looked up
- force_ccw: bool
Toggle if a counterclockwise orientation of the geometries should be enforced
- n_partitions: int
Number of workers used to lookup STARE indices in parallel
- Returns
- sids: numpy.ndarray
array of (set of) STARE index values
Examples
From points
>>> import starepandas, geopandas >>> lats = [-72.609177, -72.648590, -72.591286] >>> lons = [-41.255402, -42.054047, -41.625336] >>> geoms = geopandas.points_from_xy(lons, lats) >>> sdf = starepandas.STAREDataFrame(geometry=geoms) >>> sdf.make_sids(level=6, convex=False) 0 2299437706637111654 1 2299435211084507366 2 2299436587616075270 Name: sids, dtype: int64
From polygons
>>> gdf = geopandas.read_file(geopandas.datasets.get_path("naturalearth_lowres")) >>> sdf = starepandas.STAREDataFrame(gdf) >>> sids = sdf.make_sids(level=5)
- make_trixels(sid_column=None, n_partitions=1, wrap_lon=True, num_workers=None)#
Returns a Polygon or Multipolygon GeoSeries containing the trixels referred by the STARE indices
- Parameters
- sid_column: str
Column to use as STARE column. Default: ‘stare’
- n_partitions: int
number of (dask) workers to use to generate trixels
- wrap_lon: bool
toggle if trixels should be wraped around antimeridian.
- num_workers: int
number of workers to use
- Returns
- trixels_series: numpy.array
array of polygons or multipolygons representing the trixels
Examples
>>> import starepandas >>> sids = [648518346341351428, 900719925474099204, 1170935903116328964] >>> sdf = starepandas.STAREDataFrame(sids=sids) >>> trixels = sdf.make_trixels()
- plot(trixels=True, boundary=True, **kwargs)#
Generate a plot with matplotlib. Seminal method to GeoDataFrame.plot() All GeoDataFrame.plot() kwargs are available.
- Parameters
- trixels: bool
Toggle if trixels (rather than the SF geometry) is to be plotted
- boundary: bool
Toggle if the ring is to be plotted as a linestring rather than the polygon. Only relevant if trixels==True
Examples
>>> import starepandas >>> import geopandas >>> world = geopandas.read_file(geopandas.datasets.get_path('naturalearth_lowres')) >>> germany = world[world.name=='Germany'] >>> germany = starepandas.STAREDataFrame(germany, add_sids=True, level=8, add_trixels=True, n_partitions=1) >>> ax = germany.plot(trixels=True, boundary=True, color='y', zorder=0)
- set_sids(col, inplace=False)#
Set the StareDataFrame stare indices using either an existing column or the specified input. By default, yields a new object. The original stare column is replaced with the input.
- Parameters
- col: array-like
f stare indices or column name
- inplace: boolean
Modify the StareDataFrame in place (do not create a new object)
- Returns
- df: STAREDataFrame
the df with sids
Examples
>>> import starepandas >>> sdf = starepandas.STAREDataFrame() >>> sids = [4611686018427387903, 2299435211084507590, 2299566194809236966] >>> sdf.set_sids(sids, inplace=True)
- set_trixels(col, inplace=False)#
Set the trixel column
- Parameters
- col: array-like or string
If array like, will add the array as a new trixel column. If string, will set the df[‘col’] as the trixel column. If None, will generate trixels from the STARE column.
- inplace: bool
Modify the StareDataFrame in place (do not create a new object)
- Returns
- df: DataFrame
DataFrame or None
Examples
>>> import starepandas >>> sids = [4611686018427387903, 4611686018427387903, 4611686018427387903] >>> sdf = starepandas.STAREDataFrame(sids=sids) >>> trixels = sdf.make_trixels() >>> sdf.set_trixels(trixels, inplace=True)
- spatial_level()#
Returns the spatial level of each feature
- split_antimeridian(inplace=False, n_workers=1, drop=False)#
Splits trixels at the antimeridian
This works on trixels that cross the meridian and whose longitudes have not been wrapped around the antimeridian. I.e. when creating the trixels use sdf.make_trixels(wrap_lon=False)
Examples
>>> cities = {'name': ['midway', 'Fiji', 'Baker', 'honolulu'], ... 'lat': [28.2, -17.8, 0.2, 21.3282956], ... 'lon': [-177.35, 178.1, -176.7, -157.9]} >>> sdf = starepandas.STAREDataFrame(cities) >>> sids = starepandas.sids_from_xy(sdf.lon, sdf.lat, level=1) >>> sdf.set_sids(sids, inplace=True) >>> trixels = sdf.make_trixels(wrap_lon=False) >>> sdf.set_trixels(trixels, inplace=True) >>> cites_split = sdf.split_antimeridian(inplace=False) >>> max(max(cites_split.trixels[1].geoms[0].exterior.xy)) 180.0
- stare_disjoint(other, method='binsearch', n_workers=1)#
Returns a
Seriesofdtype('bool')with valueTruefor each geometry that is disjoint from other. This is the inverse operation of STAREDataFrame.stare_intersects()- Parameters
- other: array-like
The STARE index collection representing the spatial object to test if is intersected.
- method: str
Method for STARE intersects test ‘skiplist’, ‘binsearch’ or ‘nn’. Default: ‘binsearch’.
- n_workers: int
number of workers to be used for intersects tests
See also
STAREDataFrame.stare_intersectsintersects test
- stare_dissolve(by=None, compress_sids=True, num_workers=1, geom=False, aggfunc='first', **kwargs)#
Dissolves a dataframe subject to a field. I.e. grouping by a field/column. Seminal method to GeoDataFrame.dissolve()
- Parameters
- by: str
column to use the dissolve on. If None, dissolve all rows.
- compress_sids: bool
Toggle if STARE index values get dissolved. If not, sids will be appended. If not dissolved, there may be repetitive sids and sids that could get merged into the parent sid.
- num_workers: int
workers to use for the dissolve
- geom: bool
Toggle if the geometry column is to be dissolved. Geom column Will be dropped if set to False.
- aggfunc: str
aggregation function. E.g. ‘first’, ‘sum’, ‘mean’.
Examples
>>> import geopandas >>> world = geopandas.read_file(geopandas.datasets.get_path('naturalearth_lowres')) >>> west = world[world['continent'].isin(['Europe', 'North America'])] >>> west = starepandas.STAREDataFrame(west, add_sids=True, level=4, add_trixels=False) >>> west.stare_dissolve(by='continent', aggfunc='sum') stare ... gdp_md_est continent ... Europe [648518346341351428, 900719925474099204, 10448... ... 25284877.0 North America [1170935903116328964, 1173187702930014212, 117... ... 23505137.0
- stare_intersection(other)#
Returns a
STARESeriesof the (STARE) spatial intersection of self with other.- Parameters
- otherArray-like
The STARE index value collection representing the object to find the intersection with.
- Returns
- intersectionSTARESeries
A series of STARE index values representing the STARE interesection of each feature with other
Examples
>>> import shapely >>> nodes1 = [[102, 33], [101, 35], [105, 34], [104, 33], [102, 33]] >>> nodes2 = [[102, 34], [106, 35], [106, 33], [102, 33.5], [102, 34]] >>> polygon1 = shapely.geometry.Polygon(nodes1) >>> polygon2 = shapely.geometry.Polygon(nodes2) >>> sids1 = starepandas.sids_from_polygon(polygon1, level=5, force_ccw=True) >>> sids2 = starepandas.sids_from_polygon(polygon2, level=5, force_ccw=True)
>>> df = starepandas.STAREDataFrame(sids=[sids1]) >>> df.stare_intersection(sids2).iloc[0] array([694117292568477701, 701435641962954757, 701998591916376069])
- stare_intersects(other, method='binsearch', n_partitions=1, num_workers=None)#
Returns a
Seriesofdtype('bool')with valueTruefor each geometry that intersects other. An object is said to intersect other if its ring and interior intersects in any way with those of the other.- Parameters
- other: int or listlike
The SID collection representing the spatial object to test if is intersected.
- method: str
Method for STARE intersects test ‘skiplist’, ‘binsearch’ or ‘nn’. Default: ‘binsearch’.
- n_partitions: int
number of dask dataframe partitions to use
- num_workers: int:
number of dask workers to use
Examples
>>> germany = [4251398048237748227, 4269412446747230211, 4278419646001971203, ... 4539628424389459971, 4548635623644200963, 4566650022153682947] >>> cities = {'name': ['berlin', 'madrid'], 'sid': [4258121269174388239, 4288120002905386575]} >>> cities = starepandas.STAREDataFrame(cities, sids='sid') >>> cities.stare_intersects(germany) 0 True 1 False dtype: bool
- to_array(column, shape=None, pivot=False)#
Converts the ‘column’ to a numpy array.
Either a shape argument has to be provided or the dataframe has to contain a column x and y holding the original array coordinates.
If the dataframe has x/y columns, the column can also be pivoted. I.e. rather than reshaping according to the shape, pivoted along the x/y columns. This may be relevant if the dataframe’s row order has changed.
- Parameters
- column: str
column name to be converted to an array
- shape: tuple
x and y shape of the array. x*y has to equal the length of the dataframe
- pivot: bool
if true, rather than simple reshaping, the dataframe is pivoted along the x and y column
See also
Examples
>>> df = starepandas.STAREDataFrame({'x': [0, 0, 1, 1], ... 'y': [1, 0, 0, 1], ... 'a': [1, 2, 3, 4]}) >>> df.to_array('a', pivot=False) array([[1, 2], [3, 4]])
>>> df.to_array('a', pivot=True) array([[2, 1], [3, 4]])
- to_arrays(shape=None, pivot=False)#
Converts a STAREDataFrame into a dictionary of arrays; one array per column/field.
This may be useful to write data back to granules. Either a shape argument has to be provided or the dataframe has to contain a column x and y holding the original array coordinates. If no shape is provided, the shape is assumed to be (max(x)+1, max(y)+1).
If the dataframe has x/y columns, the column can also be pivoted. I.e. rather than reshaping according to the shape, pivoted along the x/y columns. This may be relevant if the dataframe’s row order has changed.
- Parameters
- shape: tuple
x and y shape of the array. x*y has to equal the length of the dataframe
- pivot: bool
if true, rather than simple reshaping, the dataframe is pivoted along the x and y column
See also
- to_sidecar(file_name, cover=False, shuffle=True, zlib=True)#
Writes STARE Sidecar
- to_stare_level(level, inplace=False, clear_to_level=False)#
Changes level of STARE index values to level; optionally clears location bits up to level. Caution: This method is not intended for use on feautures represented by sets of sids.
- Parameters
- inplace: bool
If True, modifies the DataFrame in place (do not create a new object).
- level: int
STARE level to change to.
- clear_to_level: bool
Toggle if the location bits below level should be cleared
- Returns
- if not inplace, returns stare index values, otherwise None
Examples
>>> sids = [2299437706637111721, 2299435211084507593, 2299566194809236969] >>> sdf = starepandas.STAREDataFrame(sids=sids) >>> sdf.to_stare_level(level=6, clear_to_level=False) sids 0 2299437706637111718 1 2299435211084507590 2 2299566194809236966
- to_stare_singllevel(level=None, inplace=False)#
Changes the STARE index values to single level representation (in contrary to multiresolution).
- Parameters
- level: int
level to change thre
- inplace: bool
If True, modifies the DataFrame in place (do not create a new object).
- Returns
- if not inplace, returns stare index values, otherwise None
Examples
>>> import geopandas >>> world = geopandas.read_file(geopandas.datasets.get_path('naturalearth_lowres')) >>> germany = world[world.name=='Germany'] >>> germany = starepandas.STAREDataFrame(germany, add_sids=True, level=6, add_trixels=False) >>> len(germany.sids.iloc[0]) 43 >>> germany_singleres = germany.to_stare_singllevel() >>> len(germany_singleres.sids.iloc[0]) 46
- trixel_area(r=None)#
Returns the approximate area of the trixel
- Parameters
- r: float or int
earth radius
- trixel_centerpoints(vertices=None)#
Returns the trixel centers as shapely points.
If vertices is set, the trixel centers are extracted from the vertices (c.f.
trixel_vertices()). If not, they are generated from the stare column.- Parameters
- vertices: tuple (vertices data structure)
If set, the centers are extracted from the vertices.
- Returns
- trixel_centerpoints: Geometery Array
Series of shapely trixel center points
Examples
>>> sids = numpy.array([4458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> centers = df.trixel_centerpoints() >>> print(centers[0]) POINT (18.4 24.09)
- trixel_centers(vertices=None)#
Returns the trixel centers.
If vertices is set, the trixel centers are extracted from the vertices (c.f.
trixel_vertices()). If not, they are generated from the stare column.- Parameters
- vertices: vertices data structure
If set, the centers are extracted from the vertices data structure.
- Returns
- trixel_centersnumpy.array
Trixel centers. First dimension are the SIDs, second dimension lon/lat.
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_centers() array([[134.9 , 80.264389]])
- trixel_centers_ecef(vertices=None)#
Returns the trixel centers as ECEF vectors.
If vertices is set, the trixel centers are extracted from the vertices (c.f.
trixel_vertices()). If not, they are generated from the stare column.- Parameters
- vertices: vertices data structure
If set, the centers are extracted from the vertices data structure.
- Returns
- trixel_centersnumpy.array
Trixel centers. First dimension are the sids, second dimension are x/y/z.
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_centers_ecef() array([[-0.11957316, 0.11957316, 0.98559856]])
- trixel_corners(vertices=None, from_trixels=False)#
Returns corners of trixels as lon/lat.
If vertices is set, the trixel corners are extracted from vertices (c.f.
trixel_vertices()). If from_trixels is True and dataframe contains trixel column, corners are extracted from trixels. If not, corners are generated from stare column- Parameters
- verticestuple (vertices data structure)
If set, the centers are extracted from the vertices.
- from_trixels: bool
If true and dataframe contains trixel column, corners are extracted from trixels.
- Returns
- cornersnumpy array
Corners of the trixels in lon/lat representation. First dimension are the SIDs, second dimension the corners (1 through 3), third dimension lon/lat.
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_corners() array([[[-170.26439001, 29.9999996 ], [ -45. , 45.00000069], [ 80.26439001, 29.9999996 ]]])
- trixel_corners_ecef(vertices=None)#
Returns ECEF norm vectors of great circles constraining the trixels.
If vertices is set, the trixel corners are extracted from vertices (c.f.
trixel_vertices()). If not, corners are generated from stare column.- Parameters
- verticestuple (vertices data structure)
If set, the centers are extracted from the vertices.
- Returns
- cornersnumpy array
Corners of the trixels in ECEF representation. First dimension are the sids, second dimension the great circles, third dimension x/y/z
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_corners_ecef() array([[[-0.85355339, -0.14644661, 0.49999999], [ 0.49999999, -0.49999999, 0.70710679], [ 0.14644661, 0.85355339, 0.49999999]]])
- trixel_grings(vertices=None)#
Returns corners of trixels as ECEF.
If vertices is set, the trixel corners are extracted from vertices (c.f.
trixel_vertices()). If not, corners are generated from stare column- Parameters
- verticestuple (vertices data structure)
If set, the centers are extracted from the vertices.
- Returns
- cornersnumpy array
ECEF norm vectors of great circles constraining the trixels. First dimension are the sids, second dimension the great circles, third dimension x/y/z
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_grings() array([[[ 0.14644661, 0.85355339, 0.49999999], [-0.85355339, -0.14644661, 0.49999999], [ 0.49999999, -0.49999999, 0.70710679]]])
- trixel_vertices()#
Returns the vertices and centerpoints of the trixels. Requires stare column to be set. Vertices are a tuple of:
the latitudes of the corners
the longitudes of the corners
the latitudes of the centers
the longitudes of the centers
- Returns
- vertices
A vertices data structure
Examples
>>> sids = numpy.array([3458764513820540928]) >>> df = starepandas.STAREDataFrame(sids=sids) >>> df.trixel_vertices() (array([29.9999996 , 45.00000069, 29.9999996 ]), array([-170.26439001, -45. , 80.26439001]), array([80.264389]), array([135.]))
- write_pods(pod_root, level, chunk_name, hex=True)#
Writes dataframe into a STAREPods hierarchy
- Parameters
- pod_root: str
Root directory of STAREPods
- level: int
level of STAREPods
- chunk_name: str
name of the pod
- hex: bool
toggle pods being hex vs int