Module phi.geom
Differentiable geometry package.
Classes:
See the phi.geom
module documentation at https://tum-pbs.github.io/PhiFlow/Geometry.html
Expand source code
"""
Differentiable geometry package.
Classes:
* `Geometry` (base type)
* `Box`
* `Sphere`
See the `phi.geom` module documentation at https://tum-pbs.github.io/PhiFlow/Geometry.html
"""
from phiml.math import stack, concat, pack_dims # for compatibility
from ._geom import Geometry, Point, assert_same_rank, invert
from ._union import union
from ._box import Box, GridCell, GridCell as UniformGrid, BaseBox, Cuboid
from ._sphere import Sphere
from ._transform import embed, infinite_cylinder
__all__ = [key for key in globals().keys() if not key.startswith('_')]
Functions
def assert_same_rank(rank1, rank2, error_message)
-
Tests that two objects have the same spatial rank. Objects can be of types:
int
,None
(no check),Geometry
,Shape
,Tensor
Expand source code
def assert_same_rank(rank1, rank2, error_message): """ Tests that two objects have the same spatial rank. Objects can be of types: `int`, `None` (no check), `Geometry`, `Shape`, `Tensor` """ rank1_, rank2_ = _rank(rank1), _rank(rank2) if rank1_ is not None and rank2_ is not None: assert rank1_ == rank2_, 'Ranks do not match: %s and %s. %s' % (rank1_, rank2_, error_message)
def concat(values: Union[tuple, list], dim: Union[str, phiml.math._shape.Shape], expand_values=False, **kwargs)
-
Concatenates a sequence of
phiml.math.magic.Shapable
objects, e.g.Tensor
, along one dimension. All values must have the same spatial, instance and channel dimensions and their sizes must be equal, except fordim
. Batch dimensions will be added as needed.Args
values
- Tuple or list of
phiml.math.magic.Shapable
, such asphiml.math.Tensor
dim
- Concatenation dimension, must be present in all
values
. The size alongdim
is determined fromvalues
and can be set to undefined (None
). expand_values
- If
True
, will first add missing dimensions to all values, not just batch dimensions. This allows tensors with different dimensions to be concatenated. The resulting tensor will have all dimensions that are present invalues
. **kwargs
- Additional keyword arguments required by specific implementations.
Adding spatial dimensions to fields requires the
bounds: Box
argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments.
Returns
Concatenated
Tensor
Examples
>>> concat([math.zeros(batch(b=10)), math.ones(batch(b=10))], 'b') (bᵇ=20) 0.500 ± 0.500 (0e+00...1e+00)
>>> concat([vec(x=1, y=0), vec(z=2.)], 'vector') (x=1.000, y=0.000, z=2.000) float64
Expand source code
def concat(values: Union[tuple, list], dim: Union[str, Shape], expand_values=False, **kwargs): """ Concatenates a sequence of `phiml.math.magic.Shapable` objects, e.g. `Tensor`, along one dimension. All values must have the same spatial, instance and channel dimensions and their sizes must be equal, except for `dim`. Batch dimensions will be added as needed. Args: values: Tuple or list of `phiml.math.magic.Shapable`, such as `phiml.math.Tensor` dim: Concatenation dimension, must be present in all `values`. The size along `dim` is determined from `values` and can be set to undefined (`None`). expand_values: If `True`, will first add missing dimensions to all values, not just batch dimensions. This allows tensors with different dimensions to be concatenated. The resulting tensor will have all dimensions that are present in `values`. **kwargs: Additional keyword arguments required by specific implementations. Adding spatial dimensions to fields requires the `bounds: Box` argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments. Returns: Concatenated `Tensor` Examples: >>> concat([math.zeros(batch(b=10)), math.ones(batch(b=10))], 'b') (bᵇ=20) 0.500 ± 0.500 (0e+00...1e+00) >>> concat([vec(x=1, y=0), vec(z=2.)], 'vector') (x=1.000, y=0.000, z=2.000) float64 """ assert len(values) > 0, f"concat() got empty sequence {values}" if isinstance(dim, Shape): dim = dim.name assert isinstance(dim, str), f"dim must be a str or Shape but got '{dim}' of type {type(dim)}" dim = auto(dim, channel).name # Add missing dimensions if expand_values: all_dims = merge_shapes(*values, allow_varying_sizes=True) all_dims = all_dims.with_dim_size(dim, 1, keep_item_names=False) values = [expand(v, all_dims.without(shape(v))) for v in values] else: for v in values: assert dim in shape(v), f"dim must be present in the shapes of all values bot got value {type(v).__name__} with shape {shape(v)}" for v in values[1:]: assert set(non_batch(v).names) == set(non_batch(values[0]).names), f"Concatenated values must have the same non-batch dimensions but got {non_batch(values[0])} and {non_batch(v)}" all_batch_dims = merge_shapes(*[shape(v).batch for v in values]) values = [expand(v, all_batch_dims) for v in values] # --- First try __concat__ --- for v in values: if isinstance(v, Shapable): if hasattr(v, '__concat__'): result = v.__concat__(values, dim, **kwargs) if result is not NotImplemented: assert isinstance(result, Shapable), f"__concat__ must return a Shapable object but got {type(result).__name__} from {type(v).__name__} {v}" return result # --- Next: try concat attributes for tree nodes --- if all(isinstance(v, PhiTreeNode) for v in values): attributes = all_attributes(values[0]) if attributes and all(all_attributes(v) == attributes for v in values): new_attrs = {} for a in attributes: common_shape = merge_shapes(*[shape(getattr(v, a)).without(dim) for v in values]) a_values = [expand(getattr(v, a), common_shape & shape(v).only(dim)) for v in values] # expand by dim if missing, and dims of others new_attrs[a] = concat(a_values, dim, **kwargs) return copy_with(values[0], **new_attrs) else: warnings.warn(f"Failed to concat values using value attributes because attributes differ among values {values}") # --- Fallback: slice and stack --- try: unstacked = sum([unstack(v, dim) for v in values], ()) except MagicNotImplemented: raise MagicNotImplemented(f"concat: No value implemented __concat__ and not all values were Sliceable along {dim}. values = {[type(v) for v in values]}") if len(unstacked) > 8: warnings.warn(f"concat() default implementation is slow on large dimensions ({dim}={len(unstacked)}). Please implement __concat__()", RuntimeWarning, stacklevel=2) dim = shape(values[0])[dim].with_size(None) try: return stack(unstacked, dim, **kwargs) except MagicNotImplemented: raise MagicNotImplemented(f"concat: No value implemented __concat__ and slices could not be stacked. values = {[type(v) for v in values]}")
def embed(geometry: phi.geom._geom.Geometry, projected_dims: Union[phiml.math._shape.Shape, str, tuple, list, None]) ‑> phi.geom._geom.Geometry
-
Adds fake spatial dimensions to a geometry. The geometry value will be constant along the added dimensions, as if it had infinite length in these directions.
Dimensions that are already present with
geometry
are ignored.Args
geometry
Geometry
projected_dims
- Additional dimensions
Returns
Geometry
with spatial rankgeometry.spatial_rank + projected_dims.rank
.Expand source code
def embed(geometry: Geometry, projected_dims: Union[math.Shape, str, tuple, list, None]) -> Geometry: """ Adds fake spatial dimensions to a geometry. The geometry value will be constant along the added dimensions, as if it had infinite length in these directions. Dimensions that are already present with `geometry` are ignored. Args: geometry: `Geometry` projected_dims: Additional dimensions Returns: `Geometry` with spatial rank `geometry.spatial_rank + projected_dims.rank`. """ if projected_dims is None: return geometry axes = parse_dim_order(projected_dims) embedded_axes = [a for a in axes if a not in geometry.shape.get_item_names('vector')] if not embedded_axes: return geometry for name in reversed(geometry.shape.get_item_names('vector')): if name not in projected_dims: axes = (name,) + axes if isinstance(geometry, BaseBox): box = geometry.corner_representation() return box * Box(**{dim: None for dim in embedded_axes}) return _EmbeddedGeometry(geometry, axes)
def infinite_cylinder(center=None, radius=None, inf_dim: Union[phiml.math._shape.Shape, tuple, list, str] = None, **center_) ‑> phi.geom._geom.Geometry
-
Creates an infinite cylinder. This is equal to embedding an
n
-dimensionalSphere
inn+1
dimensions.Args
center
- Center coordinates without
inf_dim
. Alternatively use keyword arguments. radius
- Cylinder radius.
inf_dim
- Dimension along which the cylinder is infinite.
Use
Geometry.rotated()
if the direction does not align with an axis. **center_
- Alternatively specify center coordinates without
inf_dim
as keyword arguments.
Returns
Expand source code
def infinite_cylinder(center=None, radius=None, inf_dim: Union[str, Shape, tuple, list] = None, **center_) -> Geometry: """ Creates an infinite cylinder. This is equal to embedding an `n`-dimensional `Sphere` in `n+1` dimensions. See Also: `Sphere`, `embed` Args: center: Center coordinates without `inf_dim`. Alternatively use keyword arguments. radius: Cylinder radius. inf_dim: Dimension along which the cylinder is infinite. Use `Geometry.rotated()` if the direction does not align with an axis. **center_: Alternatively specify center coordinates without `inf_dim` as keyword arguments. Returns: `Geometry` """ sphere = Sphere(center, radius, **center_) return embed(sphere, inf_dim)
def invert(geometry: phi.geom._geom.Geometry)
-
Swaps inside and outside.
Args
geometry
Geometry
to swap
Returns
New
Geometry
object with same surface but swapped normalsExpand source code
def invert(geometry: Geometry): """ Swaps inside and outside. Args: geometry: `phi.geom.Geometry` to swap Returns: New `phi.geom.Geometry` object with same surface but swapped normals """ return ~geometry
def pack_dims(value, dims: Union[str, tuple, list, set, phiml.math._shape.Shape, Callable], packed_dim: phiml.math._shape.Shape, pos: Optional[int] = None, **kwargs)
-
Compresses multiple dimensions into a single dimension by concatenating the elements. Elements along the new dimensions are laid out according to the order of
dims
. If the order ofdims
differs from the current dimension order, the tensor is transposed accordingly. This function replaces the traditionalreshape
for these cases.The type of the new dimension will be equal to the types of
dims
. Ifdims
have varying types, the new dimension will be a batch dimension.If none of
dims
exist onvalue
,packed_dim
will be added only if it is given with a definite size andvalue
is not a primitive type.See Also:
unpack_dim()
Args
value
phiml.math.magic.Shapable
, such asphiml.math.Tensor
.dims
- Dimensions to be compressed in the specified order.
packed_dim
- Single-dimension
Shape
. pos
- Index of new dimension.
None
for automatic,-1
for last,0
for first. **kwargs
- Additional keyword arguments required by specific implementations.
Adding spatial dimensions to fields requires the
bounds: Box
argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments.
Returns
Same type as
value
.Examples
>>> pack_dims(math.zeros(spatial(x=4, y=3)), spatial, instance('points')) (pointsⁱ=12) const 0.0
Expand source code
def pack_dims(value, dims: DimFilter, packed_dim: Shape, pos: Optional[int] = None, **kwargs): """ Compresses multiple dimensions into a single dimension by concatenating the elements. Elements along the new dimensions are laid out according to the order of `dims`. If the order of `dims` differs from the current dimension order, the tensor is transposed accordingly. This function replaces the traditional `reshape` for these cases. The type of the new dimension will be equal to the types of `dims`. If `dims` have varying types, the new dimension will be a batch dimension. If none of `dims` exist on `value`, `packed_dim` will be added only if it is given with a definite size and `value` is not a primitive type. See Also: `unpack_dim()` Args: value: `phiml.math.magic.Shapable`, such as `phiml.math.Tensor`. dims: Dimensions to be compressed in the specified order. packed_dim: Single-dimension `Shape`. pos: Index of new dimension. `None` for automatic, `-1` for last, `0` for first. **kwargs: Additional keyword arguments required by specific implementations. Adding spatial dimensions to fields requires the `bounds: Box` argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments. Returns: Same type as `value`. Examples: >>> pack_dims(math.zeros(spatial(x=4, y=3)), spatial, instance('points')) (pointsⁱ=12) const 0.0 """ if isinstance(value, (Number, bool)): return value assert isinstance(value, Shapable) and isinstance(value, Sliceable) and isinstance(value, Shaped), f"value must be Shapable but got {type(value)}" dims = shape(value).only(dims, reorder=True) if packed_dim in shape(value): assert packed_dim in dims, f"Cannot pack dims into new dimension {packed_dim} because it already exists on value {value} and is not packed." if len(dims) == 0 or all(dim not in shape(value) for dim in dims): return value if packed_dim.size is None else expand(value, packed_dim, **kwargs) # Inserting size=1 can cause shape errors elif len(dims) == 1: return rename_dims(value, dims, packed_dim, **kwargs) # --- First try __pack_dims__ --- if hasattr(value, '__pack_dims__'): result = value.__pack_dims__(dims.names, packed_dim, pos, **kwargs) if result is not NotImplemented: return result # --- Next try Tree Node --- if isinstance(value, PhiTreeNode): new_attributes = {a: pack_dims(getattr(value, a), dims, packed_dim, pos=pos, **kwargs) for a in all_attributes(value)} return copy_with(value, **new_attributes) # --- Fallback: unstack and stack --- if shape(value).only(dims).volume > 8: warnings.warn(f"pack_dims() default implementation is slow on large dimensions ({shape(value).only(dims)}). Please implement __pack_dims__() for {type(value).__name__} as defined in phiml.math.magic", RuntimeWarning, stacklevel=2) return stack(unstack(value, dims), packed_dim, **kwargs)
def stack(values: Union[tuple, list, dict], dim: Union[str, phiml.math._shape.Shape], expand_values=False, **kwargs)
-
Stacks
values
along the new dimensiondim
. All values must have the same spatial, instance and channel dimensions. If the dimension sizes vary, the resulting tensor will be non-uniform. Batch dimensions will be added as needed.Stacking tensors is performed lazily, i.e. the memory is allocated only when needed. This makes repeated stacking and slicing along the same dimension very efficient, i.e. jit-compiled functions will not perform these operations.
Args
values
- Collection of
phiml.math.magic.Shapable
, such asphiml.math.Tensor
If adict
, keys must be of typestr
and are used as item names alongdim
. dim
Shape
with a least one dimension. None of these dimensions can be present with any of thevalues
. Ifdim
is a single-dimension shape, its size is determined fromlen(values)
and can be left undefined (None
). Ifdim
is a multi-dimension shape, its volume must be equal tolen(values)
.expand_values
- If
True
, will first add missing dimensions to all values, not just batch dimensions. This allows tensors with different dimensions to be stacked. The resulting tensor will have all dimensions that are present invalues
. **kwargs
- Additional keyword arguments required by specific implementations.
Adding spatial dimensions to fields requires the
bounds: Box
argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments.
Returns
Tensor
containingvalues
stacked alongdim
.Examples
>>> stack({'x': 0, 'y': 1}, channel('vector')) (x=0, y=1)
>>> stack([math.zeros(batch(b=2)), math.ones(batch(b=2))], channel(c='x,y')) (x=0.000, y=1.000); (x=0.000, y=1.000) (bᵇ=2, cᶜ=x,y)
>>> stack([vec(x=1, y=0), vec(x=2, y=3.)], batch('b')) (x=1.000, y=0.000); (x=2.000, y=3.000) (bᵇ=2, vectorᶜ=x,y)
Expand source code
def stack(values: Union[tuple, list, dict], dim: Union[Shape, str], expand_values=False, **kwargs): """ Stacks `values` along the new dimension `dim`. All values must have the same spatial, instance and channel dimensions. If the dimension sizes vary, the resulting tensor will be non-uniform. Batch dimensions will be added as needed. Stacking tensors is performed lazily, i.e. the memory is allocated only when needed. This makes repeated stacking and slicing along the same dimension very efficient, i.e. jit-compiled functions will not perform these operations. Args: values: Collection of `phiml.math.magic.Shapable`, such as `phiml.math.Tensor` If a `dict`, keys must be of type `str` and are used as item names along `dim`. dim: `Shape` with a least one dimension. None of these dimensions can be present with any of the `values`. If `dim` is a single-dimension shape, its size is determined from `len(values)` and can be left undefined (`None`). If `dim` is a multi-dimension shape, its volume must be equal to `len(values)`. expand_values: If `True`, will first add missing dimensions to all values, not just batch dimensions. This allows tensors with different dimensions to be stacked. The resulting tensor will have all dimensions that are present in `values`. **kwargs: Additional keyword arguments required by specific implementations. Adding spatial dimensions to fields requires the `bounds: Box` argument specifying the physical extent of the new dimensions. Adding batch dimensions must always work without keyword arguments. Returns: `Tensor` containing `values` stacked along `dim`. Examples: >>> stack({'x': 0, 'y': 1}, channel('vector')) (x=0, y=1) >>> stack([math.zeros(batch(b=2)), math.ones(batch(b=2))], channel(c='x,y')) (x=0.000, y=1.000); (x=0.000, y=1.000) (bᵇ=2, cᶜ=x,y) >>> stack([vec(x=1, y=0), vec(x=2, y=3.)], batch('b')) (x=1.000, y=0.000); (x=2.000, y=3.000) (bᵇ=2, vectorᶜ=x,y) """ assert len(values) > 0, f"stack() got empty sequence {values}" if not dim: assert len(values) == 1, f"Only one element can be passed as `values` if no dim is passed but got {values}" from ._tensors import wrap return wrap(next(iter(values.values())) if isinstance(values, dict) else values[0]) if not isinstance(dim, Shape): dim = auto(dim) values_ = tuple(values.values()) if isinstance(values, dict) else values if not expand_values: for v in values_[1:]: assert set(non_batch(v).names) == set(non_batch(values_[0]).names), f"When expand_values=False, stacked values must have the same non-batch dimensions but got {non_batch(values_[0])} and {non_batch(v)}" # --- Add missing dimensions --- if expand_values: all_dims = merge_shapes(*values_, allow_varying_sizes=True) if isinstance(values, dict): values = {k: expand(v, all_dims.without(shape(v))) for k, v in values.items()} else: values = [expand(v, all_dims.without(shape(v))) for v in values] else: all_batch_dims = merge_shapes(*[shape(v).batch for v in values_], allow_varying_sizes=True) if isinstance(values, dict): values = {k: expand(v, all_batch_dims.without(shape(v))) for k, v in values.items()} else: values = [expand(v, all_batch_dims.without(shape(v))) for v in values] if dim.rank == 1: assert dim.size == len(values) or dim.size is None, f"stack dim size must match len(values) or be undefined but got {dim} for {len(values)} values" if dim.size is None: dim = dim.with_size(len(values)) if isinstance(values, dict): dim_item_names = tuple(values.keys()) values = tuple(values.values()) dim = dim.with_size(dim_item_names) # --- First try __stack__ --- for v in values: if hasattr(v, '__stack__'): result = v.__stack__(values, dim, **kwargs) if result is not NotImplemented: assert isinstance(result, Shapable), "__stack__ must return a Shapable object" return result # --- Next: try stacking attributes for tree nodes --- if all(isinstance(v, PhiTreeNode) for v in values): attributes = all_attributes(values[0]) if attributes and all(all_attributes(v) == attributes for v in values): new_attrs = {} for a in attributes: assert all(dim not in shape(getattr(v, a)) for v in values), f"Cannot stack attribute {a} because one values contains the stack dimension {dim}." a_values = [getattr(v, a) for v in values] if all(v is a_values[0] for v in a_values[1:]): new_attrs[a] = expand(a_values[0], dim, **kwargs) else: new_attrs[a] = stack(a_values, dim, expand_values=expand_values, **kwargs) return copy_with(values[0], **new_attrs) else: warnings.warn(f"Failed to concat values using value attributes because attributes differ among values {values}") # --- Fallback: use expand and concat --- for v in values: if not hasattr(v, '__stack__') and hasattr(v, '__concat__') and hasattr(v, '__expand__'): expanded_values = tuple([expand(v, dim.with_size(1 if dim.item_names[0] is None else dim.item_names[0][i]), **kwargs) for i, v in enumerate(values)]) if len(expanded_values) > 8: warnings.warn(f"stack() default implementation is slow on large dimensions ({dim.name}={len(expanded_values)}). Please implement __stack__()", RuntimeWarning, stacklevel=2) result = v.__concat__(expanded_values, dim.name, **kwargs) if result is not NotImplemented: assert isinstance(result, Shapable), "__concat__ must return a Shapable object" return result # --- else maybe all values are native scalars --- from ._tensors import wrap try: values = tuple([wrap(v) for v in values]) except ValueError: raise MagicNotImplemented(f"At least one item in values must be Shapable but got types {[type(v) for v in values]}") return values[0].__stack__(values, dim, **kwargs) else: # multi-dim stack assert dim.volume == len(values), f"When passing multiple stack dims, their volume must equal len(values) but got {dim} for {len(values)} values" if isinstance(values, dict): warnings.warn(f"When stacking a dict along multiple dimensions, the key names are discarded. Got keys {tuple(values.keys())}", RuntimeWarning, stacklevel=2) values = tuple(values.values()) # --- if any value implements Shapable, use stack and unpack_dim --- for v in values: if hasattr(v, '__stack__') and hasattr(v, '__unpack_dim__'): stack_dim = batch('_stack') stacked = v.__stack__(values, stack_dim, **kwargs) if stacked is not NotImplemented: assert isinstance(stacked, Shapable), "__stack__ must return a Shapable object" assert hasattr(stacked, '__unpack_dim__'), "If a value supports __unpack_dim__, the result of __stack__ must also support it." reshaped = stacked.__unpack_dim__(stack_dim.name, dim, **kwargs) if kwargs is NotImplemented: warnings.warn("__unpack_dim__ is overridden but returned NotImplemented during multi-dimensional stack. This results in unnecessary stack operations.", RuntimeWarning, stacklevel=2) else: return reshaped # --- Fallback: multi-level stack --- for dim_ in reversed(dim): values = [stack(values[i:i + dim_.size], dim_, **kwargs) for i in range(0, len(values), dim_.size)] return values[0]
def union(*geometries) ‑> phi.geom._geom.Geometry
-
Union of the given geometries. A point lies inside the union if it lies within at least one of the geometries.
Args
geometries
- arbitrary geometries with same spatial dims. Arbitrary batch dims are allowed.
*geometries
Returns
union Geometry
Expand source code
def union(*geometries) -> Geometry: """ Union of the given geometries. A point lies inside the union if it lies within at least one of the geometries. Args: geometries: arbitrary geometries with same spatial dims. Arbitrary batch dims are allowed. *geometries: Returns: union Geometry """ if len(geometries) == 1 and isinstance(geometries[0], (tuple, list)): geometries = geometries[0] if len(geometries) == 0: return NO_GEOMETRY elif len(geometries) == 1: return geometries[0] elif all(type(g) == type(geometries[0]) and isinstance(g, PhiTreeNode) for g in geometries): attrs = variable_attributes(geometries[0]) values = {a: math.stack([getattr(g, a) for g in geometries], math.instance('union')) for a in attrs} return copy_with(geometries[0], **values) else: base_geometries = () for geometry in geometries: base_geometries += geometry.geometries if isinstance(geometry, Union) else (geometry,) return Union(base_geometries)
Classes
class BaseBox
-
Abstract base type for box-like geometries.
Expand source code
class BaseBox(Geometry): # not a Subwoofer """ Abstract base type for box-like geometries. """ def __eq__(self, other): raise NotImplementedError() def __hash__(self): raise NotImplementedError() def __ne__(self, other): return not self == other @property def shape(self): raise NotImplementedError() @property def center(self) -> Tensor: raise NotImplementedError() def at(self, center: Tensor) -> 'BaseBox': return Cuboid(center, self.half_size) @property def size(self) -> Tensor: raise NotImplementedError(self) @property def half_size(self) -> Tensor: raise NotImplementedError(self) @property def lower(self) -> Tensor: raise NotImplementedError(self) @property def upper(self) -> Tensor: raise NotImplementedError(self) @property def volume(self) -> Tensor: return math.prod(self.size, 'vector') @property def shape_type(self) -> Tensor: return math.tensor('B') def bounding_radius(self): return math.vec_length(self.half_size) def bounding_half_extent(self): return self.size * 0.5 def global_to_local(self, global_position: Tensor) -> Tensor: if math.close(self.lower, 0): return global_position / self.size else: return (global_position - self.lower) / self.size def local_to_global(self, local_position): return local_position * self.size + self.lower def lies_inside(self, location): bool_inside = (location >= self.lower) & (location <= self.upper) bool_inside = math.all(bool_inside, 'vector') bool_inside = math.any(bool_inside, self.shape.instance) # union for instance dimensions return bool_inside def approximate_signed_distance(self, location: Union[Tensor, tuple]): """ Computes the signed L-infinity norm (manhattan distance) from the location to the nearest side of the box. For an outside location `l` with the closest surface point `s`, the distance is `max(abs(l - s))`. For inside locations it is `-max(abs(l - s))`. Args: location: float tensor of shape (batch_size, ..., rank) Returns: float tensor of shape (*location.shape[:-1], 1). """ center = 0.5 * (self.lower + self.upper) extent = self.upper - self.lower distance = math.abs(location - center) - extent * 0.5 distance = math.max(distance, 'vector') distance = math.min(distance, self.shape.instance) # union for instance dimensions return distance def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: loc_to_center = positions - self.center sgn_dist_from_surface = math.abs(loc_to_center) - self.half_size if outward: # --- get negative distances (particles are inside) towards the nearest boundary and add shift_amount --- distances_of_interest = (sgn_dist_from_surface == math.max(sgn_dist_from_surface, 'vector')) & (sgn_dist_from_surface < 0) shift = distances_of_interest * (sgn_dist_from_surface - shift_amount) else: shift = (sgn_dist_from_surface + shift_amount) * (sgn_dist_from_surface > 0) # get positive distances (particles are outside) and add shift_amount shift = math.where(math.abs(shift) > math.abs(loc_to_center), math.abs(loc_to_center), shift) # ensure inward shift ends at center return positions + math.where(loc_to_center < 0, 1, -1) * shift def project(self, *dimensions: str): """ Project this box into a lower-dimensional space. """ warnings.warn("Box.project(dims) is deprecated. Use Box.vector[dims] instead", DeprecationWarning, stacklevel=2) return self.vector[dimensions] def sample_uniform(self, *shape: math.Shape) -> Tensor: uniform = math.random_uniform(self.shape.non_singleton, *shape, math.channel(vector=self.spatial_rank)) return self.lower + uniform * self.size def corner_representation(self) -> 'Box': return Box(self.lower, self.upper) box = corner_representation def center_representation(self) -> 'Cuboid': return Cuboid(self.center, self.half_size) def contains(self, other: 'BaseBox'): """ Tests if the other box lies fully inside this box. """ return np.all(other.lower >= self.lower) and np.all(other.upper <= self.upper) def rotated(self, angle) -> Geometry: from ._transform import rotate return rotate(self, angle) def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return Cuboid(self.center, self.half_size * factor)
Ancestors
- phi.geom._geom.Geometry
Subclasses
- phi.geom._box.Box
- phi.geom._box.Cuboid
- phi.geom._box.GridCell
Instance variables
var center : phiml.math._tensors.Tensor
-
Center location in single channel dimension.
Expand source code
@property def center(self) -> Tensor: raise NotImplementedError()
var half_size : phiml.math._tensors.Tensor
-
Expand source code
@property def half_size(self) -> Tensor: raise NotImplementedError(self)
var lower : phiml.math._tensors.Tensor
-
Expand source code
@property def lower(self) -> Tensor: raise NotImplementedError(self)
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): raise NotImplementedError()
- A single channel dimension called
var shape_type : phiml.math._tensors.Tensor
-
Returns the type (or types) of this geometry as a string
Tensor
Boxes return'B'
, spheres return'S'
, points return'P'
. Returns'?'
for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers.Returns
String
Tensor
Expand source code
@property def shape_type(self) -> Tensor: return math.tensor('B')
var size : phiml.math._tensors.Tensor
-
Expand source code
@property def size(self) -> Tensor: raise NotImplementedError(self)
var upper : phiml.math._tensors.Tensor
-
Expand source code
@property def upper(self) -> Tensor: raise NotImplementedError(self)
var volume : phiml.math._tensors.Tensor
-
Volume of the geometry as
phiml.math.Tensor
. The result retains all batch dimensions while instance dimensions are summed over.Expand source code
@property def volume(self) -> Tensor: return math.prod(self.size, 'vector')
Methods
def approximate_signed_distance(self, location: Union[phiml.math._tensors.Tensor, tuple])
-
Computes the signed L-infinity norm (manhattan distance) from the location to the nearest side of the box. For an outside location
l
with the closest surface points
, the distance ismax(abs(l - s))
. For inside locations it is-max(abs(l - s))
.Args
location
- float tensor of shape (batch_size, …, rank)
Returns
float tensor of shape (*location.shape[:-1], 1).
Expand source code
def approximate_signed_distance(self, location: Union[Tensor, tuple]): """ Computes the signed L-infinity norm (manhattan distance) from the location to the nearest side of the box. For an outside location `l` with the closest surface point `s`, the distance is `max(abs(l - s))`. For inside locations it is `-max(abs(l - s))`. Args: location: float tensor of shape (batch_size, ..., rank) Returns: float tensor of shape (*location.shape[:-1], 1). """ center = 0.5 * (self.lower + self.upper) extent = self.upper - self.lower distance = math.abs(location - center) - extent * 0.5 distance = math.max(distance, 'vector') distance = math.min(distance, self.shape.instance) # union for instance dimensions return distance
def at(self, center: phiml.math._tensors.Tensor) ‑> phi.geom._box.BaseBox
-
Returns a copy of this
Geometry
with the center atcenter
. This is equal to callingself @ center
.See Also:
Geometry.shifted()
.Args
center
- New center as
Tensor
.
Returns
Expand source code
def at(self, center: Tensor) -> 'BaseBox': return Cuboid(center, self.half_size)
def bounding_half_extent(self)
-
The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative.
Let the bounding half-extent have value
e
in dimensiond
(extent[...,d] = e
). Then, no point of the geometry lies further away from its center point thane
alongd
(in both axis directions).:return: float vector
Args:
Returns:
Expand source code
def bounding_half_extent(self): return self.size * 0.5
def bounding_radius(self)
-
Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry.
:return: radius of type float
Args:
Returns:
Expand source code
def bounding_radius(self): return math.vec_length(self.half_size)
def box(self) ‑> phi.geom._box.Box
-
Expand source code
def corner_representation(self) -> 'Box': return Box(self.lower, self.upper)
def center_representation(self) ‑> phi.geom._box.Cuboid
-
Expand source code
def center_representation(self) -> 'Cuboid': return Cuboid(self.center, self.half_size)
def contains(self, other: BaseBox)
-
Tests if the other box lies fully inside this box.
Expand source code
def contains(self, other: 'BaseBox'): """ Tests if the other box lies fully inside this box. """ return np.all(other.lower >= self.lower) and np.all(other.upper <= self.upper)
def corner_representation(self) ‑> phi.geom._box.Box
-
Expand source code
def corner_representation(self) -> 'Box': return Box(self.lower, self.upper)
def global_to_local(self, global_position: phiml.math._tensors.Tensor) ‑> phiml.math._tensors.Tensor
-
Expand source code
def global_to_local(self, global_position: Tensor) -> Tensor: if math.close(self.lower, 0): return global_position / self.size else: return (global_position - self.lower) / self.size
def lies_inside(self, location)
-
Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside.
When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance.
Args
location
- float tensor of shape (batch_size, …, rank)
Returns
bool tensor of shape (*location.shape[:-1], 1).
Expand source code
def lies_inside(self, location): bool_inside = (location >= self.lower) & (location <= self.upper) bool_inside = math.all(bool_inside, 'vector') bool_inside = math.any(bool_inside, self.shape.instance) # union for instance dimensions return bool_inside
def local_to_global(self, local_position)
-
Expand source code
def local_to_global(self, local_position): return local_position * self.size + self.lower
def project(self, *dimensions: str)
-
Project this box into a lower-dimensional space.
Expand source code
def project(self, *dimensions: str): """ Project this box into a lower-dimensional space. """ warnings.warn("Box.project(dims) is deprecated. Use Box.vector[dims] instead", DeprecationWarning, stacklevel=2) return self.vector[dimensions]
def push(self, positions: phiml.math._tensors.Tensor, outward: bool = True, shift_amount: float = 0) ‑> phiml.math._tensors.Tensor
-
Shifts positions either into or out of geometry.
Args
positions
- Tensor holding positions to shift
outward
- Flag for indicating inward (False) or outward (True) shift
shift_amount
- Minimum distance between positions and box boundaries after shifting
Returns
Tensor holding shifted positions
Expand source code
def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: loc_to_center = positions - self.center sgn_dist_from_surface = math.abs(loc_to_center) - self.half_size if outward: # --- get negative distances (particles are inside) towards the nearest boundary and add shift_amount --- distances_of_interest = (sgn_dist_from_surface == math.max(sgn_dist_from_surface, 'vector')) & (sgn_dist_from_surface < 0) shift = distances_of_interest * (sgn_dist_from_surface - shift_amount) else: shift = (sgn_dist_from_surface + shift_amount) * (sgn_dist_from_surface > 0) # get positive distances (particles are outside) and add shift_amount shift = math.where(math.abs(shift) > math.abs(loc_to_center), math.abs(loc_to_center), shift) # ensure inward shift ends at center return positions + math.where(loc_to_center < 0, 1, -1) * shift
def rotated(self, angle) ‑> phi.geom._geom.Geometry
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle) -> Geometry: from ._transform import rotate return rotate(self, angle)
def sample_uniform(self, *shape: phiml.math._shape.Shape) ‑> phiml.math._tensors.Tensor
-
Samples uniformly distributed random points inside this volume.
Args
*shape
- How many points to sample per individual geometry.
Returns
Tensor
containing all dimensions fromGeometry.shape
,shape
as well as achannel
dimensionvector
matching the dimensionality of thisGeometry
.Expand source code
def sample_uniform(self, *shape: math.Shape) -> Tensor: uniform = math.random_uniform(self.shape.non_singleton, *shape, math.channel(vector=self.spatial_rank)) return self.lower + uniform * self.size
def scaled(self, factor: Union[phiml.math._tensors.Tensor, float]) ‑> phi.geom._geom.Geometry
-
Scales each individual geometry by
factor
. The individualcenter
points act as pivots for the operation.Args
factor: Returns:
Expand source code
def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return Cuboid(self.center, self.half_size * factor)
class Box (lower: phiml.math._tensors.Tensor = None, upper: phiml.math._tensors.Tensor = None, **size: Union[int, phiml.math._tensors.Tensor, tuple, list])
-
Simple cuboid defined by location of lower and upper corner in physical space.
Boxes can be constructed either from two positional vector arguments
(lower, upper)
or by specifying the limits by dimension name askwargs
.Examples
>>> Box(x=1, y=1) # creates a two-dimensional unit box with `lower=(0, 0)` and `upper=(1, 1)`. >>> Box(x=(None, 1), y=(0, None) # creates a Box with `lower=(-inf, 0)` and `upper=(1, inf)`.
The slicing constructor was updated in version 2.2 and now requires the dimension order as the first argument.
>>> Box['x,y', 0:1, 0:1] # creates a two-dimensional unit box with `lower=(0, 0)` and `upper=(1, 1)`. >>> Box['x,y', :1, 0:] # creates a Box with `lower=(-inf, 0)` and `upper=(1, inf)`.
Args
lower
- physical location of lower corner
upper
- physical location of upper corner
**size
- Specify size by dimension, either as
int
ortuple
containing (lower, upper).
Expand source code
class Box(BaseBox, metaclass=BoxType): """ Simple cuboid defined by location of lower and upper corner in physical space. Boxes can be constructed either from two positional vector arguments `(lower, upper)` or by specifying the limits by dimension name as `kwargs`. Examples: >>> Box(x=1, y=1) # creates a two-dimensional unit box with `lower=(0, 0)` and `upper=(1, 1)`. >>> Box(x=(None, 1), y=(0, None) # creates a Box with `lower=(-inf, 0)` and `upper=(1, inf)`. The slicing constructor was updated in version 2.2 and now requires the dimension order as the first argument. >>> Box['x,y', 0:1, 0:1] # creates a two-dimensional unit box with `lower=(0, 0)` and `upper=(1, 1)`. >>> Box['x,y', :1, 0:] # creates a Box with `lower=(-inf, 0)` and `upper=(1, inf)`. """ def __init__(self, lower: Tensor = None, upper: Tensor = None, **size: Union[int, Tensor, tuple, list]): """ Args: lower: physical location of lower corner upper: physical location of upper corner **size: Specify size by dimension, either as `int` or `tuple` containing (lower, upper). """ if lower is not None: assert isinstance(lower, Tensor), f"lower must be a Tensor but got {type(lower)}" assert 'vector' in lower.shape, "lower must have a vector dimension" assert lower.vector.item_names is not None, "vector dimension of lower must list spatial dimension order" self._lower = lower if upper is not None: assert isinstance(upper, Tensor), f"upper must be a Tensor but got {type(upper)}" assert 'vector' in upper.shape, "lower must have a vector dimension" assert upper.vector.item_names is not None, "vector dimension of lower must list spatial dimension order" self._upper = upper else: lower = [] upper = [] for item in size.values(): if isinstance(item, (tuple, list)): assert len(item) == 2, f"Box kwargs must be either dim=upper or dim=(lower,upper) but got {item}" lo, up = item lower.append(lo) upper.append(up) elif item is None: lower.append(-INF) upper.append(INF) else: lower.append(0) upper.append(item) lower = [-INF if l is None else l for l in lower] upper = [INF if u is None else u for u in upper] self._upper = math.wrap(upper, math.channel(vector=tuple(size.keys()))) self._lower = math.wrap(lower, math.channel(vector=tuple(size.keys()))) vector_shape = self._lower.shape & self._upper.shape self._lower = math.expand(self._lower, vector_shape) self._upper = math.expand(self._upper, vector_shape) if self.size.vector.item_names is None: warnings.warn("Creating a Box without item names prevents certain operations like project()", DeprecationWarning, stacklevel=2) def __getitem__(self, item): item = _keep_vector(slicing_dict(self, item)) return Box(self._lower[item], self._upper[item]) @staticmethod def __stack__(values: tuple, dim: Shape, **kwargs) -> 'Geometry': if all(isinstance(v, Box) for v in values): return NotImplemented # stack attributes else: return Geometry.__stack__(values, dim, **kwargs) def __eq__(self, other): if self._lower is None and self._upper is None: return isinstance(other, Box) return isinstance(other, BaseBox)\ and set(self.shape) == set(other.shape)\ and self.size.shape.get_size('vector') == other.size.shape.get_size('vector')\ and math.close(self._lower, other.lower)\ and math.close(self._upper, other.upper) def without(self, dims: Tuple[str, ...]): remaining = list(self.shape.get_item_names('vector')) for dim in dims: if dim in remaining: remaining.remove(dim) return self.vector[remaining] def largest(self, dim: DimFilter) -> 'Box': dim = self.shape.without('vector').only(dim) if not dim: return self return Box(math.min(self._lower, dim), math.max(self._upper, dim)) def __hash__(self): return hash(self._upper) def __variable_attrs__(self): return '_lower', '_upper' @property def shape(self): if self._lower is None or self._upper is None: return None return self._lower.shape & self._upper.shape @property def lower(self): return self._lower @property def upper(self): return self._upper @property def size(self): return self.upper - self.lower @property def center(self): return 0.5 * (self.lower + self.upper) @property def half_size(self): return self.size * 0.5 def shifted(self, delta, **delta_by_dim): return Box(self.lower + delta, self.upper + delta) def __mul__(self, other): if not isinstance(other, Box): return NotImplemented lower = self._lower.vector.unstack(self.spatial_rank) + other._lower.vector.unstack(other.spatial_rank) upper = self._upper.vector.unstack(self.spatial_rank) + other._upper.vector.unstack(other.spatial_rank) names = self._upper.vector.item_names + other._upper.vector.item_names lower = math.stack(lower, math.channel(vector=names)) upper = math.stack(upper, math.channel(vector=names)) return Box(lower, upper) def __repr__(self): if self.shape.non_channel.volume == 1: item_names = self.size.vector.item_names if item_names: return f"Box({', '.join([f'{dim}=({lo}, {up})' for dim, lo, up in zip(item_names, self._lower, self._upper)])})" else: # deprecated return 'Box[%s at %s]' % ('x'.join([str(x) for x in self.size.numpy().flatten()]), ','.join([str(x) for x in self.lower.numpy().flatten()])) else: return f'Box[shape={self.shape}]'
Ancestors
- phi.geom._box.BaseBox
- phi.geom._geom.Geometry
Instance variables
var center
-
Center location in single channel dimension.
Expand source code
@property def center(self): return 0.5 * (self.lower + self.upper)
var half_size
-
Expand source code
@property def half_size(self): return self.size * 0.5
var lower
-
Expand source code
@property def lower(self): return self._lower
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): if self._lower is None or self._upper is None: return None return self._lower.shape & self._upper.shape
- A single channel dimension called
var size
-
Expand source code
@property def size(self): return self.upper - self.lower
var upper
-
Expand source code
@property def upper(self): return self._upper
Methods
def largest(self, dim: Union[str, tuple, list, set, phiml.math._shape.Shape, Callable]) ‑> phi.geom._box.Box
-
Expand source code
def largest(self, dim: DimFilter) -> 'Box': dim = self.shape.without('vector').only(dim) if not dim: return self return Box(math.min(self._lower, dim), math.max(self._upper, dim))
def shifted(self, delta, **delta_by_dim)
-
Returns a translated version of this geometry.
See Also:
Geometry.at()
.Args
delta
- direction vector
delta
- Tensor:
Returns
Geometry
- shifted geometry
Expand source code
def shifted(self, delta, **delta_by_dim): return Box(self.lower + delta, self.upper + delta)
def without(self, dims: Tuple[str, ...])
-
Expand source code
def without(self, dims: Tuple[str, ...]): remaining = list(self.shape.get_item_names('vector')) for dim in dims: if dim in remaining: remaining.remove(dim) return self.vector[remaining]
class Cuboid (center: phiml.math._tensors.Tensor = 0, half_size: Union[phiml.math._tensors.Tensor, float] = None, **size: Union[phiml.math._tensors.Tensor, float])
-
Box specified by center position and half size.
Expand source code
class Cuboid(BaseBox): """ Box specified by center position and half size. """ def __init__(self, center: Tensor = 0, half_size: Union[float, Tensor] = None, **size: Union[float, Tensor]): if half_size is not None: assert isinstance(half_size, Tensor), "half_size must be a Tensor" assert 'vector' in half_size.shape, f"Cuboid size must have a 'vector' dimension." assert half_size.shape.get_item_names('vector') is not None, f"Vector dimension must list spatial dimensions as item names. Use the syntax Cuboid(x=x, y=y) to assign names." self._half_size = half_size else: self._half_size = math.wrap(tuple(size.values()), math.channel(vector=tuple(size.keys()))) * 0.5 center = wrap(center) if 'vector' not in center.shape or center.shape.get_item_names('vector') is None: center = math.expand(center, channel(self._half_size)) self._center = center def __eq__(self, other): if self._center is None and self._half_size is None: return isinstance(other, Cuboid) return isinstance(other, BaseBox)\ and set(self.shape) == set(other.shape)\ and math.close(self._center, other.center)\ and math.close(self._half_size, other.half_size) def __hash__(self): return hash(self._center) def __repr__(self): return f"Cuboid(center={self._center}, half_size={self._half_size})" def __getitem__(self, item): item = _keep_vector(slicing_dict(self, item)) return Cuboid(self._center[item], self._half_size[item]) @staticmethod def __stack__(values: tuple, dim: Shape, **kwargs) -> 'Geometry': if all(isinstance(v, Cuboid) for v in values): return Cuboid(math.stack([v.center for v in values], dim, **kwargs), math.stack([v.half_size for v in values], dim, **kwargs)) else: return Geometry.__stack__(values, dim, **kwargs) def __variable_attrs__(self): return '_center', '_half_size' @property def center(self): return self._center @property def half_size(self): return self._half_size @property def shape(self): if self._center is None or self._half_size is None: return None return self._center.shape & self._half_size.shape @property def size(self): return 2 * self.half_size @property def lower(self): return self.center - self.half_size @property def upper(self): return self.center + self.half_size def shifted(self, delta, **delta_by_dim) -> 'Cuboid': return Cuboid(self._center + delta, self._half_size)
Ancestors
- phi.geom._box.BaseBox
- phi.geom._geom.Geometry
Instance variables
var center
-
Center location in single channel dimension.
Expand source code
@property def center(self): return self._center
var half_size
-
Expand source code
@property def half_size(self): return self._half_size
var lower
-
Expand source code
@property def lower(self): return self.center - self.half_size
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): if self._center is None or self._half_size is None: return None return self._center.shape & self._half_size.shape
- A single channel dimension called
var size
-
Expand source code
@property def size(self): return 2 * self.half_size
var upper
-
Expand source code
@property def upper(self): return self.center + self.half_size
Methods
def shifted(self, delta, **delta_by_dim) ‑> phi.geom._box.Cuboid
-
Returns a translated version of this geometry.
See Also:
Geometry.at()
.Args
delta
- direction vector
delta
- Tensor:
Returns
Geometry
- shifted geometry
Expand source code
def shifted(self, delta, **delta_by_dim) -> 'Cuboid': return Cuboid(self._center + delta, self._half_size)
class Geometry
-
Abstract base class for N-dimensional shapes.
Main implementing classes:
- Sphere
- box family: box (generator), Box, Cuboid, BaseBox
All geometry objects support batching. Thereby any parameter defining the geometry can be varied along arbitrary batch dims. All batch dimensions are listed in Geometry.shape.
Expand source code
class Geometry: """ Abstract base class for N-dimensional shapes. Main implementing classes: * Sphere * box family: box (generator), Box, Cuboid, BaseBox All geometry objects support batching. Thereby any parameter defining the geometry can be varied along arbitrary batch dims. All batch dimensions are listed in Geometry.shape. """ @property def center(self) -> Tensor: """ Center location in single channel dimension. """ raise NotImplementedError(self) @property def shape(self) -> Shape: """ The `shape` of a `Geometry` consists of the following dimensions: * A single *channel* dimension called `'vector'` specifying the physical space * Instance dimensions denote that this geometry consists of multiple copies in the same space * Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space * Batch dimensions indicate non-interacting versions of this geometry for parallelization only. """ raise NotImplementedError() @property def volume(self) -> Tensor: """ Volume of the geometry as `phiml.math.Tensor`. The result retains all batch dimensions while instance dimensions are summed over. """ raise NotImplementedError() @property def shape_type(self) -> Tensor: """ Returns the type (or types) of this geometry as a string `Tensor` Boxes return `'B'`, spheres return `'S'`, points return `'P'`. Returns `'?'` for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers. Returns: String `Tensor` """ raise NotImplementedError() def unstack(self, dimension: str) -> tuple: """ Unstacks this Geometry along the given dimension. The shapes of the returned geometries are reduced by `dimension`. Args: dimension: dimension along which to unstack Returns: geometries: tuple of length equal to `geometry.shape.get_size(dimension)` """ return math.unstack(self, dimension) @property def spatial_rank(self) -> int: """ Number of spatial dimensions of the geometry, 1 = 1D, 2 = 2D, 3 = 3D, etc. """ return self.shape.get_size('vector') def lies_inside(self, location: Tensor) -> Tensor: """ Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside. When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance. Args: location: float tensor of shape (batch_size, ..., rank) Returns: bool tensor of shape (*location.shape[:-1], 1). """ raise NotImplementedError(self.__class__) def approximate_signed_distance(self, location: Union[Tensor, tuple]) -> Tensor: """ Computes the approximate distance from location to the surface of the geometry. Locations outside return positive values, inside negative values and zero exactly at the boundary. The exact distance metric used depends on the geometry. The approximation holds close to the surface and the distance grows to infinity as the location is moved infinitely far from the geometry. The distance metric is differentiable and its gradients are bounded at every point in space. When dealing with unions or collections of geometries (instance dimensions), the shortest distance to any instance is returned. This also holds for negative distances. Args: location: float tensor of shape (batch_size, ..., rank) location: Tensor: Returns: float tensor of shape (*location.shape[:-1], 1). """ raise NotImplementedError(self.__class__) def approximate_fraction_inside(self, other_geometry: 'Geometry', balance: Union[Tensor, Number] = 0.5) -> Tensor: """ Computes the approximate overlap between the geometry and a small other geometry. Returns 1.0 if `other_geometry` is fully enclosed in this geometry and 0.0 if there is no overlap. Close to the surface of this geometry, the fraction filled is differentiable w.r.t. the location and size of `other_geometry`. To call this method on batches of geometries of same shape, pass a batched Geometry instance. The result tensor will match the batch shape of `other_geometry`. The result may only be accurate in special cases. The given geometries may be approximated as spheres or boxes using `bounding_radius()` and `bounding_half_extent()`. The default implementation of this method approximates other_geometry as a Sphere and computes the fraction using `approximate_signed_distance()`. Args: other_geometry: `Geometry` or geometry batch for which to compute the overlap with `self`. balance: Mid-level between 0 and 1, default 0.5. This value is returned when exactly half of `other_geometry` lies inside `self`. `0.5 < balance <= 1` makes `self` seem larger while `0 <= balance < 0.5`makes `self` seem smaller. Returns: fraction of cell volume lying inside the geometry. float tensor of shape (other_geometry.batch_shape, 1). """ assert isinstance(other_geometry, Geometry) radius = other_geometry.bounding_radius() location = other_geometry.center distance = self.approximate_signed_distance(location) inside_fraction = balance - distance / radius inside_fraction = math.clip(inside_fraction, 0, 1) return inside_fraction def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: """ Shifts positions either into or out of geometry. Args: positions: Tensor holding positions to shift outward: Flag for indicating inward (False) or outward (True) shift shift_amount: Minimum distance between positions and box boundaries after shifting Returns: Tensor holding shifted positions """ raise NotImplementedError(self.__class__) def sample_uniform(self, *shape: math.Shape) -> Tensor: """ Samples uniformly distributed random points inside this volume. Args: *shape: How many points to sample per individual geometry. Returns: `Tensor` containing all dimensions from `Geometry.shape`, `shape` as well as a `channel` dimension `vector` matching the dimensionality of this `Geometry`. """ raise NotImplementedError(self.__class__) def bounding_radius(self) -> Tensor: """ Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry. :return: radius of type float Args: Returns: """ raise NotImplementedError(self.__class__) def bounding_half_extent(self) -> Tensor: """ The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative. Let the bounding half-extent have value `e` in dimension `d` (`extent[...,d] = e`). Then, no point of the geometry lies further away from its center point than `e` along `d` (in both axis directions). :return: float vector Args: Returns: """ raise NotImplementedError(self.__class__) def bounding_box(self) -> 'BaseBox': """ Returns the approximately smallest axis-aligned box that contains this `Geometry`. The center of the box may not be equal to `self.center`. Returns: `Box` or `Cuboid` that fully contains this `Geometry`. """ from ._box import Cuboid return Cuboid(self.center, half_size=self.bounding_half_extent()) def shifted(self, delta: Tensor) -> 'Geometry': """ Returns a translated version of this geometry. See Also: `Geometry.at()`. Args: delta: direction vector delta: Tensor: Returns: Geometry: shifted geometry """ return self.at(self.center + delta) def at(self, center: Tensor) -> 'Geometry': """ Returns a copy of this `Geometry` with the center at `center`. This is equal to calling `self @ center`. See Also: `Geometry.shifted()`. Args: center: New center as `Tensor`. Returns: `Geometry`. """ raise NotImplementedError def __matmul__(self, other): return self.at(other) def rotated(self, angle: Union[float, Tensor]) -> 'Geometry': """ Returns a rotated version of this geometry. The geometry is rotated about its center point. Args: angle: scalar (2d) or vector (3D+) representing delta angle Returns: Rotated `Geometry` """ raise NotImplementedError(self.__class__) def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': """ Scales each individual geometry by `factor`. The individual `center` points act as pivots for the operation. Args: factor: Returns: """ raise NotImplementedError(self.__class__) def __invert__(self): return _InvertedGeometry(self) def __eq__(self, other): """ Slow equality check. Unlike `==`, this method compares all tensor elements to check whether they are equal. Use `==` for a faster check which only checks whether the referenced tensors are the same. See Also: `shallow_equals()` """ if self is other: return True if not isinstance(other, type(self)): return False if self.shape != other.shape: return False c1 = {a: getattr(self, a) for a in variable_attributes(self)} c2 = {a: getattr(other, a) for a in variable_attributes(self)} for c in c1.keys(): if c1[c] is not c2[c] and math.any(c1[c] != c2[c]): return False return True def shallow_equals(self, other): """ Quick equality check. May return `False` even if `other == self`. However, if `True` is returned, the geometries are guaranteed to be equal. The `shallow_equals()` check does not compare all tensor elements but merely checks whether the same tensors are referenced. """ if self is other: return True if not isinstance(other, type(self)): return False if self.shape != other.shape: return False c1 = {a: getattr(self, a) for a in variable_attributes(self)} c2 = {a: getattr(other, a) for a in variable_attributes(self)} for c in c1.keys(): if c1[c] is not c2[c]: return False return True @staticmethod def __stack__(values: tuple, dim: Shape, **kwargs) -> 'Geometry': if all(type(v) == type(values[0]) for v in values): return NotImplemented # let attributes be stacked else: from ._stack import GeometryStack return GeometryStack(math.layout(values, dim)) def __flatten__(self, flat_dim: Shape, flatten_batch: bool, **kwargs) -> 'Geometry': dims = self.shape.without('vector') if not flatten_batch: dims = dims.non_batch return math.pack_dims(self, dims, flat_dim, **kwargs) def __ne__(self, other): return not self == other def __hash__(self): raise NotImplementedError(self.__class__) def __repr__(self): return f"{self.__class__.__name__}{self.shape}" def __getitem__(self, item): raise NotImplementedError # assert isinstance(item, dict), "Index must be dict of type {dim: slice/int}." # item = {dim: sel for dim, sel in item.items() if dim != 'vector'} # attrs = {a: getattr(self, a)[item] for a in variable_attributes(self)} # return copy_with(self, **attrs) def __getattr__(self, name: str) -> BoundDim: return BoundDim(self, name)
Subclasses
- phi.geom._box.BaseBox
- phi.geom._geom.Point
- phi.geom._geom._InvertedGeometry
- phi.geom._geom._NoGeometry
- phi.geom._poly_surface.PolygonSurface
- phi.geom._sphere.Sphere
- phi.geom._stack.GeometryStack
- phi.geom._transform.RotatedGeometry
- phi.geom._transform._EmbeddedGeometry
- phi.geom._union.Union
Instance variables
var center : phiml.math._tensors.Tensor
-
Center location in single channel dimension.
Expand source code
@property def center(self) -> Tensor: """ Center location in single channel dimension. """ raise NotImplementedError(self)
var shape : phiml.math._shape.Shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self) -> Shape: """ The `shape` of a `Geometry` consists of the following dimensions: * A single *channel* dimension called `'vector'` specifying the physical space * Instance dimensions denote that this geometry consists of multiple copies in the same space * Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space * Batch dimensions indicate non-interacting versions of this geometry for parallelization only. """ raise NotImplementedError()
- A single channel dimension called
var shape_type : phiml.math._tensors.Tensor
-
Returns the type (or types) of this geometry as a string
Tensor
Boxes return'B'
, spheres return'S'
, points return'P'
. Returns'?'
for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers.Returns
String
Tensor
Expand source code
@property def shape_type(self) -> Tensor: """ Returns the type (or types) of this geometry as a string `Tensor` Boxes return `'B'`, spheres return `'S'`, points return `'P'`. Returns `'?'` for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers. Returns: String `Tensor` """ raise NotImplementedError()
var spatial_rank : int
-
Number of spatial dimensions of the geometry, 1 = 1D, 2 = 2D, 3 = 3D, etc.
Expand source code
@property def spatial_rank(self) -> int: """ Number of spatial dimensions of the geometry, 1 = 1D, 2 = 2D, 3 = 3D, etc. """ return self.shape.get_size('vector')
var volume : phiml.math._tensors.Tensor
-
Volume of the geometry as
phiml.math.Tensor
. The result retains all batch dimensions while instance dimensions are summed over.Expand source code
@property def volume(self) -> Tensor: """ Volume of the geometry as `phiml.math.Tensor`. The result retains all batch dimensions while instance dimensions are summed over. """ raise NotImplementedError()
Methods
def approximate_fraction_inside(self, other_geometry: Geometry, balance: Union[phiml.math._tensors.Tensor, numbers.Number] = 0.5) ‑> phiml.math._tensors.Tensor
-
Computes the approximate overlap between the geometry and a small other geometry. Returns 1.0 if
other_geometry
is fully enclosed in this geometry and 0.0 if there is no overlap. Close to the surface of this geometry, the fraction filled is differentiable w.r.t. the location and size ofother_geometry
.To call this method on batches of geometries of same shape, pass a batched Geometry instance. The result tensor will match the batch shape of
other_geometry
.The result may only be accurate in special cases. The given geometries may be approximated as spheres or boxes using
bounding_radius()
andbounding_half_extent()
.The default implementation of this method approximates other_geometry as a Sphere and computes the fraction using
approximate_signed_distance()
.Args
other_geometry
Geometry
or geometry batch for which to compute the overlap withself
.balance
- Mid-level between 0 and 1, default 0.5.
This value is returned when exactly half of
other_geometry
lies insideself
.0.5 < balance <= 1
makesself
seem larger while0 <= balance < 0.5
makesself
seem smaller.
Returns
fraction of cell volume lying inside the geometry. float tensor of shape (other_geometry.batch_shape, 1).
Expand source code
def approximate_fraction_inside(self, other_geometry: 'Geometry', balance: Union[Tensor, Number] = 0.5) -> Tensor: """ Computes the approximate overlap between the geometry and a small other geometry. Returns 1.0 if `other_geometry` is fully enclosed in this geometry and 0.0 if there is no overlap. Close to the surface of this geometry, the fraction filled is differentiable w.r.t. the location and size of `other_geometry`. To call this method on batches of geometries of same shape, pass a batched Geometry instance. The result tensor will match the batch shape of `other_geometry`. The result may only be accurate in special cases. The given geometries may be approximated as spheres or boxes using `bounding_radius()` and `bounding_half_extent()`. The default implementation of this method approximates other_geometry as a Sphere and computes the fraction using `approximate_signed_distance()`. Args: other_geometry: `Geometry` or geometry batch for which to compute the overlap with `self`. balance: Mid-level between 0 and 1, default 0.5. This value is returned when exactly half of `other_geometry` lies inside `self`. `0.5 < balance <= 1` makes `self` seem larger while `0 <= balance < 0.5`makes `self` seem smaller. Returns: fraction of cell volume lying inside the geometry. float tensor of shape (other_geometry.batch_shape, 1). """ assert isinstance(other_geometry, Geometry) radius = other_geometry.bounding_radius() location = other_geometry.center distance = self.approximate_signed_distance(location) inside_fraction = balance - distance / radius inside_fraction = math.clip(inside_fraction, 0, 1) return inside_fraction
def approximate_signed_distance(self, location: Union[phiml.math._tensors.Tensor, tuple]) ‑> phiml.math._tensors.Tensor
-
Computes the approximate distance from location to the surface of the geometry. Locations outside return positive values, inside negative values and zero exactly at the boundary.
The exact distance metric used depends on the geometry. The approximation holds close to the surface and the distance grows to infinity as the location is moved infinitely far from the geometry. The distance metric is differentiable and its gradients are bounded at every point in space.
When dealing with unions or collections of geometries (instance dimensions), the shortest distance to any instance is returned. This also holds for negative distances.
Args
location
- float tensor of shape (batch_size, …, rank)
location
- Tensor:
Returns
float tensor of shape (*location.shape[:-1], 1).
Expand source code
def approximate_signed_distance(self, location: Union[Tensor, tuple]) -> Tensor: """ Computes the approximate distance from location to the surface of the geometry. Locations outside return positive values, inside negative values and zero exactly at the boundary. The exact distance metric used depends on the geometry. The approximation holds close to the surface and the distance grows to infinity as the location is moved infinitely far from the geometry. The distance metric is differentiable and its gradients are bounded at every point in space. When dealing with unions or collections of geometries (instance dimensions), the shortest distance to any instance is returned. This also holds for negative distances. Args: location: float tensor of shape (batch_size, ..., rank) location: Tensor: Returns: float tensor of shape (*location.shape[:-1], 1). """ raise NotImplementedError(self.__class__)
def at(self, center: phiml.math._tensors.Tensor) ‑> phi.geom._geom.Geometry
-
Returns a copy of this
Geometry
with the center atcenter
. This is equal to callingself @ center
.See Also:
Geometry.shifted()
.Args
center
- New center as
Tensor
.
Returns
Expand source code
def at(self, center: Tensor) -> 'Geometry': """ Returns a copy of this `Geometry` with the center at `center`. This is equal to calling `self @ center`. See Also: `Geometry.shifted()`. Args: center: New center as `Tensor`. Returns: `Geometry`. """ raise NotImplementedError
def bounding_box(self) ‑> BaseBox
-
Returns the approximately smallest axis-aligned box that contains this
Geometry
. The center of the box may not be equal toself.center
.Returns
Expand source code
def bounding_box(self) -> 'BaseBox': """ Returns the approximately smallest axis-aligned box that contains this `Geometry`. The center of the box may not be equal to `self.center`. Returns: `Box` or `Cuboid` that fully contains this `Geometry`. """ from ._box import Cuboid return Cuboid(self.center, half_size=self.bounding_half_extent())
def bounding_half_extent(self) ‑> phiml.math._tensors.Tensor
-
The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative.
Let the bounding half-extent have value
e
in dimensiond
(extent[...,d] = e
). Then, no point of the geometry lies further away from its center point thane
alongd
(in both axis directions).:return: float vector
Args:
Returns:
Expand source code
def bounding_half_extent(self) -> Tensor: """ The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative. Let the bounding half-extent have value `e` in dimension `d` (`extent[...,d] = e`). Then, no point of the geometry lies further away from its center point than `e` along `d` (in both axis directions). :return: float vector Args: Returns: """ raise NotImplementedError(self.__class__)
def bounding_radius(self) ‑> phiml.math._tensors.Tensor
-
Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry.
:return: radius of type float
Args:
Returns:
Expand source code
def bounding_radius(self) -> Tensor: """ Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry. :return: radius of type float Args: Returns: """ raise NotImplementedError(self.__class__)
def lies_inside(self, location: phiml.math._tensors.Tensor) ‑> phiml.math._tensors.Tensor
-
Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside.
When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance.
Args
location
- float tensor of shape (batch_size, …, rank)
Returns
bool tensor of shape (*location.shape[:-1], 1).
Expand source code
def lies_inside(self, location: Tensor) -> Tensor: """ Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside. When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance. Args: location: float tensor of shape (batch_size, ..., rank) Returns: bool tensor of shape (*location.shape[:-1], 1). """ raise NotImplementedError(self.__class__)
def push(self, positions: phiml.math._tensors.Tensor, outward: bool = True, shift_amount: float = 0) ‑> phiml.math._tensors.Tensor
-
Shifts positions either into or out of geometry.
Args
positions
- Tensor holding positions to shift
outward
- Flag for indicating inward (False) or outward (True) shift
shift_amount
- Minimum distance between positions and box boundaries after shifting
Returns
Tensor holding shifted positions
Expand source code
def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: """ Shifts positions either into or out of geometry. Args: positions: Tensor holding positions to shift outward: Flag for indicating inward (False) or outward (True) shift shift_amount: Minimum distance between positions and box boundaries after shifting Returns: Tensor holding shifted positions """ raise NotImplementedError(self.__class__)
def rotated(self, angle: Union[phiml.math._tensors.Tensor, float]) ‑> phi.geom._geom.Geometry
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle: Union[float, Tensor]) -> 'Geometry': """ Returns a rotated version of this geometry. The geometry is rotated about its center point. Args: angle: scalar (2d) or vector (3D+) representing delta angle Returns: Rotated `Geometry` """ raise NotImplementedError(self.__class__)
def sample_uniform(self, *shape: phiml.math._shape.Shape) ‑> phiml.math._tensors.Tensor
-
Samples uniformly distributed random points inside this volume.
Args
*shape
- How many points to sample per individual geometry.
Returns
Tensor
containing all dimensions fromGeometry.shape
,shape
as well as achannel
dimensionvector
matching the dimensionality of thisGeometry
.Expand source code
def sample_uniform(self, *shape: math.Shape) -> Tensor: """ Samples uniformly distributed random points inside this volume. Args: *shape: How many points to sample per individual geometry. Returns: `Tensor` containing all dimensions from `Geometry.shape`, `shape` as well as a `channel` dimension `vector` matching the dimensionality of this `Geometry`. """ raise NotImplementedError(self.__class__)
def scaled(self, factor: Union[phiml.math._tensors.Tensor, float]) ‑> phi.geom._geom.Geometry
-
Scales each individual geometry by
factor
. The individualcenter
points act as pivots for the operation.Args
factor: Returns:
Expand source code
def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': """ Scales each individual geometry by `factor`. The individual `center` points act as pivots for the operation. Args: factor: Returns: """ raise NotImplementedError(self.__class__)
def shallow_equals(self, other)
-
Quick equality check. May return
False
even ifother == self
. However, ifTrue
is returned, the geometries are guaranteed to be equal.The
shallow_equals()
check does not compare all tensor elements but merely checks whether the same tensors are referenced.Expand source code
def shallow_equals(self, other): """ Quick equality check. May return `False` even if `other == self`. However, if `True` is returned, the geometries are guaranteed to be equal. The `shallow_equals()` check does not compare all tensor elements but merely checks whether the same tensors are referenced. """ if self is other: return True if not isinstance(other, type(self)): return False if self.shape != other.shape: return False c1 = {a: getattr(self, a) for a in variable_attributes(self)} c2 = {a: getattr(other, a) for a in variable_attributes(self)} for c in c1.keys(): if c1[c] is not c2[c]: return False return True
def shifted(self, delta: phiml.math._tensors.Tensor) ‑> phi.geom._geom.Geometry
-
Returns a translated version of this geometry.
See Also:
Geometry.at()
.Args
delta
- direction vector
delta
- Tensor:
Returns
Geometry
- shifted geometry
Expand source code
def shifted(self, delta: Tensor) -> 'Geometry': """ Returns a translated version of this geometry. See Also: `Geometry.at()`. Args: delta: direction vector delta: Tensor: Returns: Geometry: shifted geometry """ return self.at(self.center + delta)
def unstack(self, dimension: str) ‑> tuple
-
Unstacks this Geometry along the given dimension. The shapes of the returned geometries are reduced by
dimension
.Args
dimension
- dimension along which to unstack
Returns
geometries
- tuple of length equal to
geometry.shape.get_size(dimension)
Expand source code
def unstack(self, dimension: str) -> tuple: """ Unstacks this Geometry along the given dimension. The shapes of the returned geometries are reduced by `dimension`. Args: dimension: dimension along which to unstack Returns: geometries: tuple of length equal to `geometry.shape.get_size(dimension)` """ return math.unstack(self, dimension)
class GridCell (resolution: phiml.math._shape.Shape, bounds: phi.geom._box.BaseBox)
-
An instance of GridCell represents all cells of a regular grid as a batch of boxes.
Expand source code
class GridCell(BaseBox): """ An instance of GridCell represents all cells of a regular grid as a batch of boxes. """ def __init__(self, resolution: math.Shape, bounds: BaseBox): assert resolution.spatial_rank == resolution.rank, f"resolution must be purely spatial but got {resolution}" assert resolution.spatial_rank == bounds.spatial_rank, f"bounds must match dimensions of resolution but got {bounds} for resolution {resolution}" assert set(bounds.vector.item_names) == set(resolution.names) self._resolution = resolution.only(bounds.vector.item_names, reorder=True) self._bounds = bounds self._shape = self._resolution & bounds.shape.non_spatial @property def resolution(self): return self._resolution @property def bounds(self): return self._bounds @property def spatial_rank(self) -> int: return self._resolution.spatial_rank @property def center(self): local_coords = math.meshgrid(**{dim.name: math.linspace(0.5 / dim.size, 1 - 0.5 / dim.size, dim) for dim in self.resolution}) points = self.bounds.local_to_global(local_coords) return points @property def grid_size(self): return self._bounds.size @property def size(self): return self.bounds.size / math.wrap(self.resolution.sizes) @property def dx(self): return self.bounds.size / self.resolution @property def lower(self): return self.center - self.half_size @property def upper(self): return self.center + self.half_size @property def half_size(self): return self.bounds.size / self.resolution.sizes / 2 def __getitem__(self, item): item = slicing_dict(self, item) bounds = self._bounds dx = self.size gather_dict = {} for dim, selection in item.items(): if dim in self._resolution: if isinstance(selection, int): start = selection stop = selection + 1 elif isinstance(selection, slice): start = selection.start or 0 if start < 0: start += self.resolution.get_size(dim) stop = selection.stop or self.resolution.get_size(dim) if stop < 0: stop += self.resolution.get_size(dim) assert selection.step is None or selection.step == 1 else: raise ValueError(f"Illegal selection: {item}") dim_mask = math.wrap(self.resolution.mask(dim)) lower = bounds.lower + start * dim_mask * dx upper = bounds.upper + (stop - self.resolution.get_size(dim)) * dim_mask * dx bounds = Box(lower, upper) gather_dict[dim] = slice(start, stop) resolution = self._resolution.after_gather(gather_dict) return GridCell(resolution, bounds[{d: s for d, s in item.items() if d != 'vector'}]) def __pack_dims__(self, dims: Tuple[str, ...], packed_dim: Shape, pos: Union[int, None], **kwargs) -> 'Cuboid': return math.pack_dims(self.center_representation(), dims, packed_dim, pos, **kwargs) @staticmethod def __stack__(values: tuple, dim: Shape, **kwargs) -> 'Geometry': from ._stack import GeometryStack return GeometryStack(math.layout(values, dim)) def list_cells(self, dim_name): center = math.pack_dims(self.center, self._shape.spatial.names, dim_name) return Cuboid(center, self.half_size) def stagger(self, dim: str, lower: bool, upper: bool): dim_mask = np.array(self.resolution.mask(dim)) unit = self.bounds.size / self.resolution * dim_mask bounds = Box(self.bounds.lower + unit * (-0.5 if lower else 0.5), self.bounds.upper + unit * (0.5 if upper else -0.5)) ext_res = self.resolution.sizes + dim_mask * (int(lower) + int(upper) - 1) return GridCell(self.resolution.with_sizes(ext_res), bounds) def padded(self, widths: dict): resolution, bounds = self.resolution, self.bounds for dim, (lower, upper) in widths.items(): masked_dx = self.dx * math.dim_mask(self.resolution, dim) resolution = resolution.with_dim_size(dim, self.resolution.get_size(dim) + lower + upper) bounds = Box(bounds.lower - masked_dx * lower, bounds.upper + masked_dx * upper) return GridCell(resolution, bounds) # def face_centers(self, staggered_name='staggered'): # face_centers = [self.extend_symmetric(dim).center for dim in self.shape.spatial.names] # return math.channel_stack(face_centers, staggered_name) @property def shape(self): return self._shape def shifted(self, delta: Tensor, **delta_by_dim) -> BaseBox: # delta += math.padded_stack() if delta.shape.spatial_rank == 0: return GridCell(self.resolution, self.bounds.shifted(delta)) else: center = self.center + delta return Cuboid(center, self.half_size) def rotated(self, angle) -> Geometry: raise NotImplementedError("Grids cannot be rotated. Use center_representation() to convert it to Cuboids first.") def __eq__(self, other): return isinstance(other, GridCell) and self._bounds == other._bounds and self._resolution == other._resolution def shallow_equals(self, other): return self == other def __hash__(self): return hash(self._resolution) + hash(self._bounds) def __repr__(self): return f"{self._resolution}, bounds={self._bounds}" def __variable_attrs__(self): return '_center', '_half_size' def __with_attrs__(self, **attrs): return copy_with(self.center_representation(), **attrs) @property def _center(self): return self.center @property def _half_size(self): return self.half_size
Ancestors
- phi.geom._box.BaseBox
- phi.geom._geom.Geometry
Instance variables
var bounds
-
Expand source code
@property def bounds(self): return self._bounds
var center
-
Center location in single channel dimension.
Expand source code
@property def center(self): local_coords = math.meshgrid(**{dim.name: math.linspace(0.5 / dim.size, 1 - 0.5 / dim.size, dim) for dim in self.resolution}) points = self.bounds.local_to_global(local_coords) return points
var dx
-
Expand source code
@property def dx(self): return self.bounds.size / self.resolution
var grid_size
-
Expand source code
@property def grid_size(self): return self._bounds.size
var half_size
-
Expand source code
@property def half_size(self): return self.bounds.size / self.resolution.sizes / 2
var lower
-
Expand source code
@property def lower(self): return self.center - self.half_size
var resolution
-
Expand source code
@property def resolution(self): return self._resolution
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): return self._shape
- A single channel dimension called
var size
-
Expand source code
@property def size(self): return self.bounds.size / math.wrap(self.resolution.sizes)
var spatial_rank : int
-
Number of spatial dimensions of the geometry, 1 = 1D, 2 = 2D, 3 = 3D, etc.
Expand source code
@property def spatial_rank(self) -> int: return self._resolution.spatial_rank
var upper
-
Expand source code
@property def upper(self): return self.center + self.half_size
Methods
def list_cells(self, dim_name)
-
Expand source code
def list_cells(self, dim_name): center = math.pack_dims(self.center, self._shape.spatial.names, dim_name) return Cuboid(center, self.half_size)
def padded(self, widths: dict)
-
Expand source code
def padded(self, widths: dict): resolution, bounds = self.resolution, self.bounds for dim, (lower, upper) in widths.items(): masked_dx = self.dx * math.dim_mask(self.resolution, dim) resolution = resolution.with_dim_size(dim, self.resolution.get_size(dim) + lower + upper) bounds = Box(bounds.lower - masked_dx * lower, bounds.upper + masked_dx * upper) return GridCell(resolution, bounds)
def rotated(self, angle) ‑> phi.geom._geom.Geometry
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle) -> Geometry: raise NotImplementedError("Grids cannot be rotated. Use center_representation() to convert it to Cuboids first.")
def shallow_equals(self, other)
-
Quick equality check. May return
False
even ifother == self
. However, ifTrue
is returned, the geometries are guaranteed to be equal.The
shallow_equals()
check does not compare all tensor elements but merely checks whether the same tensors are referenced.Expand source code
def shallow_equals(self, other): return self == other
def shifted(self, delta: phiml.math._tensors.Tensor, **delta_by_dim) ‑> phi.geom._box.BaseBox
-
Returns a translated version of this geometry.
See Also:
Geometry.at()
.Args
delta
- direction vector
delta
- Tensor:
Returns
Geometry
- shifted geometry
Expand source code
def shifted(self, delta: Tensor, **delta_by_dim) -> BaseBox: # delta += math.padded_stack() if delta.shape.spatial_rank == 0: return GridCell(self.resolution, self.bounds.shifted(delta)) else: center = self.center + delta return Cuboid(center, self.half_size)
def stagger(self, dim: str, lower: bool, upper: bool)
-
Expand source code
def stagger(self, dim: str, lower: bool, upper: bool): dim_mask = np.array(self.resolution.mask(dim)) unit = self.bounds.size / self.resolution * dim_mask bounds = Box(self.bounds.lower + unit * (-0.5 if lower else 0.5), self.bounds.upper + unit * (0.5 if upper else -0.5)) ext_res = self.resolution.sizes + dim_mask * (int(lower) + int(upper) - 1) return GridCell(self.resolution.with_sizes(ext_res), bounds)
class UniformGrid (resolution: phiml.math._shape.Shape, bounds: phi.geom._box.BaseBox)
-
An instance of GridCell represents all cells of a regular grid as a batch of boxes.
Expand source code
class GridCell(BaseBox): """ An instance of GridCell represents all cells of a regular grid as a batch of boxes. """ def __init__(self, resolution: math.Shape, bounds: BaseBox): assert resolution.spatial_rank == resolution.rank, f"resolution must be purely spatial but got {resolution}" assert resolution.spatial_rank == bounds.spatial_rank, f"bounds must match dimensions of resolution but got {bounds} for resolution {resolution}" assert set(bounds.vector.item_names) == set(resolution.names) self._resolution = resolution.only(bounds.vector.item_names, reorder=True) self._bounds = bounds self._shape = self._resolution & bounds.shape.non_spatial @property def resolution(self): return self._resolution @property def bounds(self): return self._bounds @property def spatial_rank(self) -> int: return self._resolution.spatial_rank @property def center(self): local_coords = math.meshgrid(**{dim.name: math.linspace(0.5 / dim.size, 1 - 0.5 / dim.size, dim) for dim in self.resolution}) points = self.bounds.local_to_global(local_coords) return points @property def grid_size(self): return self._bounds.size @property def size(self): return self.bounds.size / math.wrap(self.resolution.sizes) @property def dx(self): return self.bounds.size / self.resolution @property def lower(self): return self.center - self.half_size @property def upper(self): return self.center + self.half_size @property def half_size(self): return self.bounds.size / self.resolution.sizes / 2 def __getitem__(self, item): item = slicing_dict(self, item) bounds = self._bounds dx = self.size gather_dict = {} for dim, selection in item.items(): if dim in self._resolution: if isinstance(selection, int): start = selection stop = selection + 1 elif isinstance(selection, slice): start = selection.start or 0 if start < 0: start += self.resolution.get_size(dim) stop = selection.stop or self.resolution.get_size(dim) if stop < 0: stop += self.resolution.get_size(dim) assert selection.step is None or selection.step == 1 else: raise ValueError(f"Illegal selection: {item}") dim_mask = math.wrap(self.resolution.mask(dim)) lower = bounds.lower + start * dim_mask * dx upper = bounds.upper + (stop - self.resolution.get_size(dim)) * dim_mask * dx bounds = Box(lower, upper) gather_dict[dim] = slice(start, stop) resolution = self._resolution.after_gather(gather_dict) return GridCell(resolution, bounds[{d: s for d, s in item.items() if d != 'vector'}]) def __pack_dims__(self, dims: Tuple[str, ...], packed_dim: Shape, pos: Union[int, None], **kwargs) -> 'Cuboid': return math.pack_dims(self.center_representation(), dims, packed_dim, pos, **kwargs) @staticmethod def __stack__(values: tuple, dim: Shape, **kwargs) -> 'Geometry': from ._stack import GeometryStack return GeometryStack(math.layout(values, dim)) def list_cells(self, dim_name): center = math.pack_dims(self.center, self._shape.spatial.names, dim_name) return Cuboid(center, self.half_size) def stagger(self, dim: str, lower: bool, upper: bool): dim_mask = np.array(self.resolution.mask(dim)) unit = self.bounds.size / self.resolution * dim_mask bounds = Box(self.bounds.lower + unit * (-0.5 if lower else 0.5), self.bounds.upper + unit * (0.5 if upper else -0.5)) ext_res = self.resolution.sizes + dim_mask * (int(lower) + int(upper) - 1) return GridCell(self.resolution.with_sizes(ext_res), bounds) def padded(self, widths: dict): resolution, bounds = self.resolution, self.bounds for dim, (lower, upper) in widths.items(): masked_dx = self.dx * math.dim_mask(self.resolution, dim) resolution = resolution.with_dim_size(dim, self.resolution.get_size(dim) + lower + upper) bounds = Box(bounds.lower - masked_dx * lower, bounds.upper + masked_dx * upper) return GridCell(resolution, bounds) # def face_centers(self, staggered_name='staggered'): # face_centers = [self.extend_symmetric(dim).center for dim in self.shape.spatial.names] # return math.channel_stack(face_centers, staggered_name) @property def shape(self): return self._shape def shifted(self, delta: Tensor, **delta_by_dim) -> BaseBox: # delta += math.padded_stack() if delta.shape.spatial_rank == 0: return GridCell(self.resolution, self.bounds.shifted(delta)) else: center = self.center + delta return Cuboid(center, self.half_size) def rotated(self, angle) -> Geometry: raise NotImplementedError("Grids cannot be rotated. Use center_representation() to convert it to Cuboids first.") def __eq__(self, other): return isinstance(other, GridCell) and self._bounds == other._bounds and self._resolution == other._resolution def shallow_equals(self, other): return self == other def __hash__(self): return hash(self._resolution) + hash(self._bounds) def __repr__(self): return f"{self._resolution}, bounds={self._bounds}" def __variable_attrs__(self): return '_center', '_half_size' def __with_attrs__(self, **attrs): return copy_with(self.center_representation(), **attrs) @property def _center(self): return self.center @property def _half_size(self): return self.half_size
Ancestors
- phi.geom._box.BaseBox
- phi.geom._geom.Geometry
Instance variables
var bounds
-
Expand source code
@property def bounds(self): return self._bounds
var center
-
Center location in single channel dimension.
Expand source code
@property def center(self): local_coords = math.meshgrid(**{dim.name: math.linspace(0.5 / dim.size, 1 - 0.5 / dim.size, dim) for dim in self.resolution}) points = self.bounds.local_to_global(local_coords) return points
var dx
-
Expand source code
@property def dx(self): return self.bounds.size / self.resolution
var grid_size
-
Expand source code
@property def grid_size(self): return self._bounds.size
var half_size
-
Expand source code
@property def half_size(self): return self.bounds.size / self.resolution.sizes / 2
var lower
-
Expand source code
@property def lower(self): return self.center - self.half_size
var resolution
-
Expand source code
@property def resolution(self): return self._resolution
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): return self._shape
- A single channel dimension called
var size
-
Expand source code
@property def size(self): return self.bounds.size / math.wrap(self.resolution.sizes)
var spatial_rank : int
-
Number of spatial dimensions of the geometry, 1 = 1D, 2 = 2D, 3 = 3D, etc.
Expand source code
@property def spatial_rank(self) -> int: return self._resolution.spatial_rank
var upper
-
Expand source code
@property def upper(self): return self.center + self.half_size
Methods
def list_cells(self, dim_name)
-
Expand source code
def list_cells(self, dim_name): center = math.pack_dims(self.center, self._shape.spatial.names, dim_name) return Cuboid(center, self.half_size)
def padded(self, widths: dict)
-
Expand source code
def padded(self, widths: dict): resolution, bounds = self.resolution, self.bounds for dim, (lower, upper) in widths.items(): masked_dx = self.dx * math.dim_mask(self.resolution, dim) resolution = resolution.with_dim_size(dim, self.resolution.get_size(dim) + lower + upper) bounds = Box(bounds.lower - masked_dx * lower, bounds.upper + masked_dx * upper) return GridCell(resolution, bounds)
def rotated(self, angle) ‑> phi.geom._geom.Geometry
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle) -> Geometry: raise NotImplementedError("Grids cannot be rotated. Use center_representation() to convert it to Cuboids first.")
def shallow_equals(self, other)
-
Quick equality check. May return
False
even ifother == self
. However, ifTrue
is returned, the geometries are guaranteed to be equal.The
shallow_equals()
check does not compare all tensor elements but merely checks whether the same tensors are referenced.Expand source code
def shallow_equals(self, other): return self == other
def shifted(self, delta: phiml.math._tensors.Tensor, **delta_by_dim) ‑> phi.geom._box.BaseBox
-
Returns a translated version of this geometry.
See Also:
Geometry.at()
.Args
delta
- direction vector
delta
- Tensor:
Returns
Geometry
- shifted geometry
Expand source code
def shifted(self, delta: Tensor, **delta_by_dim) -> BaseBox: # delta += math.padded_stack() if delta.shape.spatial_rank == 0: return GridCell(self.resolution, self.bounds.shifted(delta)) else: center = self.center + delta return Cuboid(center, self.half_size)
def stagger(self, dim: str, lower: bool, upper: bool)
-
Expand source code
def stagger(self, dim: str, lower: bool, upper: bool): dim_mask = np.array(self.resolution.mask(dim)) unit = self.bounds.size / self.resolution * dim_mask bounds = Box(self.bounds.lower + unit * (-0.5 if lower else 0.5), self.bounds.upper + unit * (0.5 if upper else -0.5)) ext_res = self.resolution.sizes + dim_mask * (int(lower) + int(upper) - 1) return GridCell(self.resolution.with_sizes(ext_res), bounds)
class Point (location: phiml.math._tensors.Tensor)
-
Points have zero volume and are determined by a single location. An instance of
Point
represents a single n-dimensional point or a batch of points.Expand source code
class Point(Geometry): """ Points have zero volume and are determined by a single location. An instance of `Point` represents a single n-dimensional point or a batch of points. """ def __init__(self, location: math.Tensor): assert 'vector' in location.shape, "location must have a vector dimension" assert location.shape.get_item_names('vector') is not None, "Vector dimension needs to list spatial dimension as item names." self._location = location @property def center(self) -> Tensor: return self._location @property def shape(self) -> Shape: return self._location.shape def unstack(self, dimension: str) -> tuple: return tuple(Point(loc) for loc in self._location.unstack(dimension)) def lies_inside(self, location: Tensor) -> Tensor: return expand(math.wrap(False), shape(location).without('vector')) def approximate_signed_distance(self, location: Union[Tensor, tuple]) -> Tensor: return math.vec_abs(location - self._location) def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: return positions def bounding_radius(self) -> Tensor: return math.zeros() def bounding_half_extent(self) -> Tensor: return math.zeros() def at(self, center: Tensor) -> 'Geometry': return Point(center) def rotated(self, angle) -> 'Geometry': return self def __hash__(self): return hash(self._location) def __variable_attrs__(self): return '_location', @property def volume(self) -> Tensor: return math.wrap(0) @property def shape_type(self) -> Tensor: return math.tensor('P') def sample_uniform(self, *shape: math.Shape) -> Tensor: raise NotImplementedError def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return self def __getitem__(self, item): return Point(self._location[_keep_vector(slicing_dict(self, item))])
Ancestors
- phi.geom._geom.Geometry
Instance variables
var center : phiml.math._tensors.Tensor
-
Center location in single channel dimension.
Expand source code
@property def center(self) -> Tensor: return self._location
var shape : phiml.math._shape.Shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self) -> Shape: return self._location.shape
- A single channel dimension called
var shape_type : phiml.math._tensors.Tensor
-
Returns the type (or types) of this geometry as a string
Tensor
Boxes return'B'
, spheres return'S'
, points return'P'
. Returns'?'
for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers.Returns
String
Tensor
Expand source code
@property def shape_type(self) -> Tensor: return math.tensor('P')
var volume : phiml.math._tensors.Tensor
-
Volume of the geometry as
phiml.math.Tensor
. The result retains all batch dimensions while instance dimensions are summed over.Expand source code
@property def volume(self) -> Tensor: return math.wrap(0)
Methods
def approximate_signed_distance(self, location: Union[phiml.math._tensors.Tensor, tuple]) ‑> phiml.math._tensors.Tensor
-
Computes the approximate distance from location to the surface of the geometry. Locations outside return positive values, inside negative values and zero exactly at the boundary.
The exact distance metric used depends on the geometry. The approximation holds close to the surface and the distance grows to infinity as the location is moved infinitely far from the geometry. The distance metric is differentiable and its gradients are bounded at every point in space.
When dealing with unions or collections of geometries (instance dimensions), the shortest distance to any instance is returned. This also holds for negative distances.
Args
location
- float tensor of shape (batch_size, …, rank)
location
- Tensor:
Returns
float tensor of shape (*location.shape[:-1], 1).
Expand source code
def approximate_signed_distance(self, location: Union[Tensor, tuple]) -> Tensor: return math.vec_abs(location - self._location)
def at(self, center: phiml.math._tensors.Tensor) ‑> phi.geom._geom.Geometry
-
Returns a copy of this
Geometry
with the center atcenter
. This is equal to callingself @ center
.See Also:
Geometry.shifted()
.Args
center
- New center as
Tensor
.
Returns
Expand source code
def at(self, center: Tensor) -> 'Geometry': return Point(center)
def bounding_half_extent(self) ‑> phiml.math._tensors.Tensor
-
The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative.
Let the bounding half-extent have value
e
in dimensiond
(extent[...,d] = e
). Then, no point of the geometry lies further away from its center point thane
alongd
(in both axis directions).:return: float vector
Args:
Returns:
Expand source code
def bounding_half_extent(self) -> Tensor: return math.zeros()
def bounding_radius(self) ‑> phiml.math._tensors.Tensor
-
Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry.
:return: radius of type float
Args:
Returns:
Expand source code
def bounding_radius(self) -> Tensor: return math.zeros()
def lies_inside(self, location: phiml.math._tensors.Tensor) ‑> phiml.math._tensors.Tensor
-
Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside.
When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance.
Args
location
- float tensor of shape (batch_size, …, rank)
Returns
bool tensor of shape (*location.shape[:-1], 1).
Expand source code
def lies_inside(self, location: Tensor) -> Tensor: return expand(math.wrap(False), shape(location).without('vector'))
def push(self, positions: phiml.math._tensors.Tensor, outward: bool = True, shift_amount: float = 0) ‑> phiml.math._tensors.Tensor
-
Shifts positions either into or out of geometry.
Args
positions
- Tensor holding positions to shift
outward
- Flag for indicating inward (False) or outward (True) shift
shift_amount
- Minimum distance between positions and box boundaries after shifting
Returns
Tensor holding shifted positions
Expand source code
def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: return positions
def rotated(self, angle) ‑> phi.geom._geom.Geometry
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle) -> 'Geometry': return self
def sample_uniform(self, *shape: phiml.math._shape.Shape) ‑> phiml.math._tensors.Tensor
-
Samples uniformly distributed random points inside this volume.
Args
*shape
- How many points to sample per individual geometry.
Returns
Tensor
containing all dimensions fromGeometry.shape
,shape
as well as achannel
dimensionvector
matching the dimensionality of thisGeometry
.Expand source code
def sample_uniform(self, *shape: math.Shape) -> Tensor: raise NotImplementedError
def scaled(self, factor: Union[phiml.math._tensors.Tensor, float]) ‑> phi.geom._geom.Geometry
-
Scales each individual geometry by
factor
. The individualcenter
points act as pivots for the operation.Args
factor: Returns:
Expand source code
def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return self
def unstack(self, dimension: str) ‑> tuple
-
Unstacks this Geometry along the given dimension. The shapes of the returned geometries are reduced by
dimension
.Args
dimension
- dimension along which to unstack
Returns
geometries
- tuple of length equal to
geometry.shape.get_size(dimension)
Expand source code
def unstack(self, dimension: str) -> tuple: return tuple(Point(loc) for loc in self._location.unstack(dimension))
class Sphere (center: phiml.math._tensors.Tensor = None, radius: Union[phiml.math._tensors.Tensor, float] = None, **center_: Union[phiml.math._tensors.Tensor, float])
-
N-dimensional sphere. Defined through center position and radius.
Args
center
- Sphere center as
Tensor
withvector
dimension. The spatial dimension order should be specified in thevector
dimension via item names. radius
- Sphere radius as
float
orTensor
**center_
- Specifies center when the
center
argument is not given. Center position by dimension, e.g.x=0.5, y=0.2
.
Expand source code
class Sphere(Geometry): """ N-dimensional sphere. Defined through center position and radius. """ def __init__(self, center: Tensor = None, radius: Union[float, Tensor] = None, **center_: Union[float, Tensor]): """ Args: center: Sphere center as `Tensor` with `vector` dimension. The spatial dimension order should be specified in the `vector` dimension via item names. radius: Sphere radius as `float` or `Tensor` **center_: Specifies center when the `center` argument is not given. Center position by dimension, e.g. `x=0.5, y=0.2`. """ if center is not None: assert isinstance(center, Tensor), "center must be a Tensor" assert 'vector' in center.shape, f"Sphere center must have a 'vector' dimension." assert center.shape.get_item_names('vector') is not None, f"Vector dimension must list spatial dimensions as item names. Use the syntax Sphere(x=x, y=y) to assign names." self._center = center else: self._center = wrap(tuple(center_.values()), math.channel(vector=tuple(center_.keys()))) assert radius is not None, "radius must be specified." self._radius = wrap(radius) assert 'vector' not in self._radius.shape, f"Sphere radius must not vary along vector but got {radius}" @property def shape(self): if self._center is None or self._radius is None: return None return self._center.shape & self._radius.shape @property def radius(self): return self._radius @property def center(self): return self._center @property def volume(self) -> math.Tensor: if self.spatial_rank == 1: return 2 * self._radius elif self.spatial_rank == 2: return math.PI * self._radius ** 2 elif self.spatial_rank == 3: return 4 / 3 * math.PI * self._radius ** 3 else: raise NotImplementedError() # n = self.spatial_rank # return math.pi ** (n // 2) / math.faculty(math.ceil(n / 2)) * self._radius ** n @property def shape_type(self) -> Tensor: return math.tensor('S') def lies_inside(self, location): distance_squared = math.sum((location - self.center) ** 2, dim='vector') return math.any(distance_squared <= self.radius ** 2, self.shape.instance) # union for instance dimensions def approximate_signed_distance(self, location: Union[Tensor, tuple]): """ Computes the exact distance from location to the closest point on the sphere. Very close to the sphere center, the distance takes a constant value. Args: location: float tensor of shape (batch_size, ..., rank) Returns: float tensor of shape (*location.shape[:-1], 1). """ distance_squared = math.vec_squared(location - self.center) distance_squared = math.maximum(distance_squared, self.radius * 1e-2) # Prevent infinite spatial_gradient at sphere center distance = math.sqrt(distance_squared) return math.min(distance - self.radius, self.shape.instance) # union for instance dimensions def sample_uniform(self, *shape: math.Shape): raise NotImplementedError('Not yet implemented') # ToDo def bounding_radius(self): return self.radius def bounding_half_extent(self): return expand(self.radius, self._center.shape.only('vector')) def at(self, center: Tensor) -> 'Geometry': return Sphere(center, self._radius) def rotated(self, angle): return self def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return Sphere(self.center, self.radius * factor) def __variable_attrs__(self): return '_center', '_radius' def __getitem__(self, item): item = slicing_dict(self, item) return Sphere(self._center[_keep_vector(item)], self._radius[item]) def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: raise NotImplementedError() def __hash__(self): return hash(self._center) + hash(self._radius)
Ancestors
- phi.geom._geom.Geometry
Instance variables
var center
-
Center location in single channel dimension.
Expand source code
@property def center(self): return self._center
var radius
-
Expand source code
@property def radius(self): return self._radius
var shape
-
The
shape
of aGeometry
consists of the following dimensions:- A single channel dimension called
'vector'
specifying the physical space - Instance dimensions denote that this geometry consists of multiple copies in the same space
- Spatial dimensions denote a crystal (repeating structure) of this geometric primitive in space
- Batch dimensions indicate non-interacting versions of this geometry for parallelization only.
Expand source code
@property def shape(self): if self._center is None or self._radius is None: return None return self._center.shape & self._radius.shape
- A single channel dimension called
var shape_type : phiml.math._tensors.Tensor
-
Returns the type (or types) of this geometry as a string
Tensor
Boxes return'B'
, spheres return'S'
, points return'P'
. Returns'?'
for unknown types, e.g. a union over multiple types. Custom types can return their own identifiers.Returns
String
Tensor
Expand source code
@property def shape_type(self) -> Tensor: return math.tensor('S')
var volume : phiml.math._tensors.Tensor
-
Volume of the geometry as
phiml.math.Tensor
. The result retains all batch dimensions while instance dimensions are summed over.Expand source code
@property def volume(self) -> math.Tensor: if self.spatial_rank == 1: return 2 * self._radius elif self.spatial_rank == 2: return math.PI * self._radius ** 2 elif self.spatial_rank == 3: return 4 / 3 * math.PI * self._radius ** 3 else: raise NotImplementedError() # n = self.spatial_rank # return math.pi ** (n // 2) / math.faculty(math.ceil(n / 2)) * self._radius ** n
Methods
def approximate_signed_distance(self, location: Union[phiml.math._tensors.Tensor, tuple])
-
Computes the exact distance from location to the closest point on the sphere. Very close to the sphere center, the distance takes a constant value.
Args
location
- float tensor of shape (batch_size, …, rank)
Returns
float tensor of shape (*location.shape[:-1], 1).
Expand source code
def approximate_signed_distance(self, location: Union[Tensor, tuple]): """ Computes the exact distance from location to the closest point on the sphere. Very close to the sphere center, the distance takes a constant value. Args: location: float tensor of shape (batch_size, ..., rank) Returns: float tensor of shape (*location.shape[:-1], 1). """ distance_squared = math.vec_squared(location - self.center) distance_squared = math.maximum(distance_squared, self.radius * 1e-2) # Prevent infinite spatial_gradient at sphere center distance = math.sqrt(distance_squared) return math.min(distance - self.radius, self.shape.instance) # union for instance dimensions
def at(self, center: phiml.math._tensors.Tensor) ‑> phi.geom._geom.Geometry
-
Returns a copy of this
Geometry
with the center atcenter
. This is equal to callingself @ center
.See Also:
Geometry.shifted()
.Args
center
- New center as
Tensor
.
Returns
Expand source code
def at(self, center: Tensor) -> 'Geometry': return Sphere(center, self._radius)
def bounding_half_extent(self)
-
The bounding half-extent sets a limit on the outer-most point for each coordinate axis. Each component is non-negative.
Let the bounding half-extent have value
e
in dimensiond
(extent[...,d] = e
). Then, no point of the geometry lies further away from its center point thane
alongd
(in both axis directions).:return: float vector
Args:
Returns:
Expand source code
def bounding_half_extent(self): return expand(self.radius, self._center.shape.only('vector'))
def bounding_radius(self)
-
Returns the radius of a Sphere object that fully encloses this geometry. The sphere is centered at the center of this geometry.
:return: radius of type float
Args:
Returns:
Expand source code
def bounding_radius(self): return self.radius
def lies_inside(self, location)
-
Tests whether the given location lies inside or outside of the geometry. Locations on the surface count as inside.
When dealing with unions or collections of geometries (instance dimensions), a point lies inside the geometry if it lies inside any instance.
Args
location
- float tensor of shape (batch_size, …, rank)
Returns
bool tensor of shape (*location.shape[:-1], 1).
Expand source code
def lies_inside(self, location): distance_squared = math.sum((location - self.center) ** 2, dim='vector') return math.any(distance_squared <= self.radius ** 2, self.shape.instance) # union for instance dimensions
def push(self, positions: phiml.math._tensors.Tensor, outward: bool = True, shift_amount: float = 0) ‑> phiml.math._tensors.Tensor
-
Shifts positions either into or out of geometry.
Args
positions
- Tensor holding positions to shift
outward
- Flag for indicating inward (False) or outward (True) shift
shift_amount
- Minimum distance between positions and box boundaries after shifting
Returns
Tensor holding shifted positions
Expand source code
def push(self, positions: Tensor, outward: bool = True, shift_amount: float = 0) -> Tensor: raise NotImplementedError()
def rotated(self, angle)
-
Returns a rotated version of this geometry. The geometry is rotated about its center point.
Args
angle
- scalar (2d) or vector (3D+) representing delta angle
Returns
Rotated
Geometry
Expand source code
def rotated(self, angle): return self
def sample_uniform(self, *shape: phiml.math._shape.Shape)
-
Samples uniformly distributed random points inside this volume.
Args
*shape
- How many points to sample per individual geometry.
Returns
Tensor
containing all dimensions fromGeometry.shape
,shape
as well as achannel
dimensionvector
matching the dimensionality of thisGeometry
.Expand source code
def sample_uniform(self, *shape: math.Shape): raise NotImplementedError('Not yet implemented') # ToDo
def scaled(self, factor: Union[phiml.math._tensors.Tensor, float]) ‑> phi.geom._geom.Geometry
-
Scales each individual geometry by
factor
. The individualcenter
points act as pivots for the operation.Args
factor: Returns:
Expand source code
def scaled(self, factor: Union[float, Tensor]) -> 'Geometry': return Sphere(self.center, self.radius * factor)