Modulus Sym Geometry
Defines base class for all geometries
- class modulus.sym.geometry.geometry.Geometry(curves, sdf, dims, bounds, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>, interior_epsilon=1e-06)[source]
Bases:
object
Base class for all geometries
- property dims
returns: dims – output can be [‘x’], [‘x’,’y’], or [‘x’,’y’,’z’] :rtype: List[srt]
- repeat(spacing: float, repeat_lower: List[int], repeat_higher: List[int], center: Union[None, List[float]] = None)[source]
Finite Repetition of geometry.
- Parameters
spacing (float) – Spacing between each repetition.
repeat_lower (List[int]) – How many repetitions going in negative direction.
repeat_upper (List[int]) – How many repetitions going in positive direction.
center (Union[None, List[Union[float, sympy.Basic]]] = None) – If given then center the rotation around this point.
-
rotate(angle: ~typing.Union[float, ~sympy.core.basic.Basic], axis: str = 'z', center: ~typing.Union[None, ~typing.List[float]] = None, parameterization=
)[source] Rotates geometry.
- Parameters
angle (Union[float, sympy.Basic]) – Angle of rotate in radians. Can be a sympy expression if parameterizing.
axis (str) – Axis of rotation. Default is “z”.
center (Union[None, List[Union[float, sympy.Basic]]] = None) – If given then center the rotation around this point.
parameterization (Parameterization) – Parameterization if translation is parameterized.
- sample_boundary(nr_points: int, criteria: Optional[Basic] = None, parameterization: Optional[Parameterization] = None, quasirandom: bool = False)[source]
Samples the surface or perimeter of the geometry.
- Parameters
nr_points (int) – number of points to sample on boundary.
criteria (Union[sympy.Basic, None]) – Only sample points that satisfy this criteria.
parameterization (Union[Parameterization, None], optional) – If the geometry is parameterized then you can provide ranges for the parameters with this. By default the sampling will be done with the internal parameterization.
quasirandom (bool) – If true then sample the points using the Halton sequences. Default is False.
- Returns
- Return type
points – Dictionary contain a point cloud sampled uniformly. For example in 2D it would be ``` points = {‘x’: np.ndarray (N, 1),
’y’: np.ndarray (N, 1), ‘normal_x’: np.ndarray (N, 1), ‘normal_y’: np.ndarray (N, 1), ‘area’: np.ndarray (N, 1)}
``` The area value can be used for Monte Carlo integration like the following, total_area = np.sum(points[‘area’])
Dict[str, np.ndarray]
- sample_interior(nr_points: int, bounds: Optional[Bounds] = None, criteria: Optional[Basic] = None, parameterization: Optional[Parameterization] = None, compute_sdf_derivatives: bool = False, quasirandom: bool = False)[source]
Samples the interior of the geometry.
- Parameters
nr_points (int) – number of points to sample.
bounds (Union[Bounds, None]) – Bounds to sample points from. For example, bounds = Bounds({Parameter(‘x’): (0, 1), Parameter(‘y’): (0, 1)}). By default the internal bounds will be used.
criteria (Union[sympy.Basic, None]) – Only sample points that satisfy this criteria.
parameterization (Union[Parameterization, None]) – If the geometry is parameterized then you can provide ranges for the parameters with this.
compute_sdf_derivatives (bool) – Compute sdf derivatives if true.
quasirandom (bool) – If true then sample the points using the Halton sequences. Default is False.
- Returns
- Return type
points – Dictionary contain a point cloud sampled uniformly. For example in 2D it would be ``` points = {‘x’: np.ndarray (N, 1),
’y’: np.ndarray (N, 1), ‘sdf’: np.ndarray (N, 1), ‘area’: np.ndarray (N, 1)}
``` The area value can be used for Monte Carlo integration like the following, total_area = np.sum(points[‘area’])
Dict[str, np.ndarray]
-
scale(x: ~typing.Union[float, ~sympy.core.basic.Basic], parameterization: ~modulus.sym.geometry.parameterization.Parameterization =
)[source] Scales geometry.
- Parameters
x (Union[float, sympy.Basic]) – Scale factor. Can be a sympy expression if parameterizing.
parameterization (Parameterization) – Parameterization if scale factor is parameterized.
-
translate(xyz: ~typing.List[~typing.Union[float, ~sympy.core.basic.Basic]], parameterization: ~modulus.sym.geometry.parameterization.Parameterization =
)[source] Translates geometry.
- Parameters
xyz (List[Union[float, sympy.Basic]]) – Translation. Can be a sympy expression if parameterizing.
parameterization (Parameterization) – Parameterization if translation is parameterized.
- class modulus.sym.geometry.adf.ADF[source]
Bases:
Module
Used for hard imposition of boundary conditions. Currently supports 2d geometries and Dirichlet boundary conditions. Contributors: M. A. Nabian, R. Gladstone, H. Meidani, N. Sukumar, A. Srivastava Reference: “Sukumar, N. and Srivastava, A., 2021.
Exact imposition of boundary conditions with distance functions in physics-informed deep neural networks. Computer Methods in Applied Mechanics and Engineering, p.114333.”
- static circle_adf(points: Tuple[Tensor], radius: float, center: Tuple[float]) → Tensor[source]
Computes the pointwise approximate distance for a circle
- Parameters
points (Tuple[torch.Tensor]) – ADF will be computed on these points
radius (float) – Radius of the circle
center (Tuple[float]) – Center of the circle
- Returns
- Return type
omega – pointwise approximate distance
torch.Tensor
- forward(invar)[source]
Defines the computation performed at every call.
Should be overridden by all subclasses.
- static infinite_line_adf(points: Tuple[Tensor], point_1: Tuple[float], point_2: Tuple[float]) → Tensor[source]
Computes the pointwise approximate distance for an infinite line
- Parameters
points (Tuple[torch.Tensor]) – ADF will be computed on these points
point_1 (Tuple[float]) – One of the two points that form the infinite line
point_2 (Tuple[float]) – One of the two points that form the infinite line
- Returns
- Return type
omega – pointwise approximate distance
torch.Tensor
- static line_segment_adf(points: Tuple[Tensor], point_1: Tuple[float], point_2: Tuple[float]) → Tensor[source]
Computes the pointwise approximate distance for a line segment
- Parameters
points (Tuple[torch.Tensor]) – ADF will be computed on these points
point_1 (Tuple[float]) – Point on one end of the line segment
point_2 (Tuple[float]) – Point on the other ned of the line segment
- Returns
- Return type
omega – pointwise approximate distance
torch.Tensor
- static r_equivalence(omegas: List[Tensor], m: float = 2.0) → Tensor[source]
Computes the R-equivalence of a collection of approximate distance functions
- Parameters
omegas (List[torch.Tensor]) – List of ADFs used to compute the R-equivalence.
m (float) – Normalization order
- Returns
- Return type
omega_E – R-equivalence distance
torch.Tensor
- static transfinite_interpolation(bases: List[Tensor], indx: int, eps: float = 1e-08) → Tensor[source]
Performs transfinite interpolation of the boundary conditions
- Parameters
bases (List[torch.Tensor]) – List of ADFs used for the transfinite interpolation.
indx (int) – index of the interpolation basis
eps (float) – Small value to avoid division by zero
- Returns
- Return type
w – Interpolation basis corresponding to the input index
torch.Tensor
- static trimmed_circle_adf(points: Tuple[Tensor], point_1: Tuple[float], point_2: Tuple[float], sign: int, radius: float, center: float) → Tensor[source]
Computes the pointwise approximate distance of a trimmed circle
- Parameters
points (Tuple[torch.Tensor]) – ADF will be computed on these points
point_1 (Tuple[float]) – One of the two points that form the trimming infinite line
point_2 (Tuple[float]) – One of the two points that form the trimming infinite line
sign (int) – Specifies the trimming side
radius (float) – Radius of the circle
center (Tuple[float]) – Center of the circle
- Returns
- Return type
omega – pointwise approximate distance
torch.Tensor
Defines different Curve objects
- class modulus.sym.geometry.curve.Curve(sample, dims, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
object
A Curve object that keeps track of the surface/perimeter of a geometry. The curve object also contains normals and area/length of curve.
- approx_area(parameterization=<modulus.sym.geometry.parameterization.Parameterization object>, criteria=None, approx_nr=10000, quasirandom=False)[source]
- Parameters
parameterization (dict with of Parameters and their ranges) – If the curve is parameterized then you can provide ranges for the parameters with this.
criteria (None, SymPy boolean exprs) – Calculate area discarding regions that don’t satisfy this criteria.
approx_nr (int) – Area might be difficult to compute if parameterized. In this case the area is approximated by sampling area, approx_nr number of times. This amounts to monte carlo integration.
- Returns
- Return type
area – area of curve
float
- property dims
returns: dims – output can be [‘x’], [‘x’,’y’], or [‘x’,’y’,’z’] :rtype: list of strings
- rotate(angle, axis, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
rotate curve
- Parameters
x (float, SymPy Symbol/Exprs) – scale factor.
- scale(x, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
scale curve
- Parameters
x (float, SymPy Symbol/Exprs) – scale factor.
- translate(xyz, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
translate curve
- Parameters
xyz (tuple of floats, ints, SymPy Symbol/Exprs) – translate curve by these values.
- class modulus.sym.geometry.curve.SympyCurve(functions, parameterization, area, criteria=None)[source]
Bases:
Curve
Curve defined by sympy functions
- Parameters
functions (dictionary of SymPy Exprs) –
Parameterized curve in 1, 2 or 3 dimensions. For example, a circle might have:
TODO: refactor to remove normals.
ranges (dictionary of Sympy Symbols and ranges) – This gives the ranges for the parameters in the parameterized curve. For example, a circle might have ranges = {theta: (0, 2*pi)}.
area (float, int, SymPy Exprs) – The surface area/perimeter of the curve.
criteria (SymPy Boolean Function) – If this boolean expression is false then we do not sample their on curve. This can be used to enforce uniform sample probability.
-
class modulus.sym.geometry.parameterization.Bounds(bound_ranges: ~typing.Dict[~modulus.sym.geometry.parameterization.Parameter, ~typing.Tuple[~typing.Union[float, ~sympy.core.basic.Basic], ~typing.Union[float, ~sympy.core.basic.Basic]]], parameterization: ~modulus.sym.geometry.parameterization.Parameterization =
)[source] Bases:
object
A object used to store bounds for geometries.
- Parameters
bound_ranges (Dict[Parameter, Tuple[Union[float, sympy.Basic], Union[float, sympy.Basic]]) – Dictionary of Parameters with names “x”, “y”, or “z”. The value given for each of these is a tuple of the lower and upper bound. Sympy expressions can be used to define these upper and lower bounds.
parameterization (Parameterization) – A Parameterization object used when the upper and lower bounds are parameterized.
- property dims
returns: dims – output can be [‘x’], [‘x’,’y’], or [‘x’,’y’,’z’] :rtype: list of strings
- sample(nr_points: int, parameterization: Optional[Parameterization] = None, quasirandom: bool = False)[source]
Sample points in Bounds.
- Parameters
nr_points (int) – Number of points sampled from parameterization.
parameterization (Parameterization) – Given if sampling bounds with different parameterization then the internal one stored in Bounds. Default is to not use this.
quasirandom (bool) – If true then sample the points using Halton sequences. Default is False.
- volume(parameterization: Optional[Parameterization] = None)[source]
Compute volume of bounds.
- Parameters
parameterization (Parameterization) – Given if sampling bounds with different parameterization then the internal one stored in Bounds. Default is to not use this.
- class modulus.sym.geometry.parameterization.OrderedParameterization(param_ranges, key)[source]
Bases:
Parameterization
A object used to store ordered parameterization information about user-specified keys.
- Parameters
param_ranges (Dict[Parameter, Union[float, Tuple[float, float], np.ndarray (N, 1)]) – Dictionary of Parameters and their ranges. The ranges can be one of the following types, :obj: Float will sample the parameter equal to this value. :obj: Tuple of two float as the bounding range to sample parameter from. :obj: np.ndarray as a discrete list of possible values for the parameter.
- sample(nr_points: int, quasirandom: bool = False, sort: Optional = 'ascending')[source]
Sample ordered parameterization values.
- Parameters
nr_points (int) – Number of points sampled from parameterization.
quasirandom (bool) – If true then sample the points using Halton sequences. Default is False.
sort (None or {'ascending','descending'}) – If ‘ascending’ then sample the sorted points in ascending order. If ‘descending’ then sample the sorted points in descending order. Default is ‘ascending’.
- class modulus.sym.geometry.parameterization.Parameter(name: str)[source]
Bases:
Symbol
A Symbolic object used to parameterize geometries.
Currently this only overloads the Sympy Symbol class however capabilities may be expanded in the future.
- Parameters
name (str) – Name given to parameter.
- class modulus.sym.geometry.parameterization.Parameterization(param_ranges: Dict[Parameter, Union[float, Tuple[float, float], ndarray]] = {})[source]
Bases:
object
A object used to store parameterization information about geometries.
- Parameters
param_ranges (Dict[Parameter, Union[float, Tuple[float, float], np.ndarray (N, 1)]) – Dictionary of Parameters and their ranges. The ranges can be one of the following types, :obj: Float will sample the parameter equal to this value. :obj: Tuple of two float as the bounding range to sample parameter from. :obj: np.ndarray as a discrete list of possible values for the parameter.
- sample(nr_points: int, quasirandom: bool = False)[source]
Sample parameterization values.
- Parameters
nr_points (int) – Number of points sampled from parameterization.
quasirandom (bool) – If true then sample the points using Halton sequences. Default is False.
- class modulus.sym.geometry.primitives_1d.Line1D(point_1, point_2, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
1D Line along x-axis
- Parameters
point_1 (int or float) – lower bound point of line
point_2 (int or float) – upper bound point of line
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_1d.Point1D(point, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
1D Point along x-axis
- Parameters
point (int or float) – x coordinate of the point
parameterization (Parameterization) – Parameterization of geometry.
Primitives for 2D geometries see https://www.iquilezles.org/www/articles/distfunctions/distfunctions.html
- class modulus.sym.geometry.primitives_2d.Channel2D(point_1, point_2, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Channel (no bounding curves in x-direction)
- Parameters
point_1 (tuple with 2 ints or floats) – lower bound point of channel
point_2 (tuple with 2 ints or floats) – upper bound point of channel
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Circle(center, radius, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Circle
- Parameters
center (tuple with 2 ints or floats) – center point of circle
radius (int or float) – radius of circle
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Ellipse(center, major, minor, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Ellipse
- Parameters
center (tuple with 2 ints or floats) – center point of circle
radius (int or float) – radius of circle
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Line(point_1, point_2, normal=1, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Line parallel to y-axis
- Parameters
point_1 (tuple with 2 ints or floats) – lower bound point of line segment
point_2 (tuple with 2 ints or floats) – upper bound point of line segment
normal (int or float) – normal direction of line (+1 or -1)
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Polygon(points, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Polygon
- Parameters
points (list of tuple with 2 ints or floats) – lower bound point of line segment
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Rectangle(point_1, point_2, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Rectangle
- Parameters
point_1 (tuple with 2 ints or floats) – lower bound point of rectangle
point_2 (tuple with 2 ints or floats) – upper bound point of rectangle
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_2d.Triangle(center, base, height, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Isosceles Triangle Symmetrical axis parallel to y-axis
- Parameters
center (tuple with 2 ints or floats) – center of base of triangle
base (int or float) – base of triangle
height (int or float) – height of triangle
parameterization (Parameterization) – Parameterization of geometry.
Primitives for 3D geometries see https://www.iquilezles.org/www/articles/distfunctions/distfunctions.html
- class modulus.sym.geometry.primitives_3d.Box(point_1, point_2, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Box/Cuboid
- Parameters
point_1 (tuple with 3 ints or floats) – lower bound point of box
point_2 (tuple with 3 ints or floats) – upper bound point of box
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Channel(point_1, point_2, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Channel (no bounding surfaces in x-direction)
- Parameters
point_1 (tuple with 3 ints or floats) – lower bound point of channel
point_2 (tuple with 3 ints or floats) – upper bound point of channel
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Cone(center, radius, height, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Cone Axis parallel to z-axis
- Parameters
center (tuple with 3 ints or floats) – base center of cone
radius (int or float) – base radius of cone
height (int or float) – height of cone
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Cylinder(center, radius, height, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Cylinder Axis parallel to z-axis
- Parameters
center (tuple with 3 ints or floats) – center of cylinder
radius (int or float) – radius of cylinder
height (int or float) – height of cylinder
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.ElliCylinder(center, a, b, height, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Elliptical Cylinder Axis parallel to z-axis
Approximation based on 4-arc ellipse construction https://www.researchgate.net/publication/241719740_Approximating_an_ellipse_with_four_circular_arcs
Please manually ensure a>b
- Parameters
center (tuple with 3 ints or floats) – center of base of ellipse
a (int or float) – semi-major axis of ellipse
b (int or float) – semi-minor axis of ellipse
height (int or float) – height of elliptical cylinder
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.IsoTriangularPrism(center, base, height, height_prism, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
2D Isosceles Triangular Prism Symmetrical axis parallel to y-axis
- Parameters
center (tuple with 3 ints or floats) – center of base of triangle
base (int or float) – base of triangle
height (int or float) – height of triangle
height_prism (int or float) – height of triangular prism
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Plane(point_1, point_2, normal=1, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Plane perpendicular to x-axis
- Parameters
point_1 (tuple with 3 ints or floats) – lower bound point of plane
point_2 (tuple with 3 ints or floats) – upper bound point of plane
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Sphere(center, radius, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Sphere
- Parameters
center (tuple with 3 ints or floats) – center of sphere
radius (int or float) – radius of sphere
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Tetrahedron(center, radius, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Tetrahedron The 4 symmetrically placed points are on a unit sphere. Centroid of the tetrahedron is at origin and lower face is parallel to x-y plane Reference: https://en.wikipedia.org/wiki/Tetrahedron
- Parameters
center (tuple with 3 ints or floats) – centroid of tetrahedron
radius (int or float) – radius of circumscribed sphere
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.Torus(center, radius, radius_tube, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Torus
- Parameters
center (tuple with 3 ints or floats) – center of torus
radius (int or float) – distance from center to center of tube (major radius)
radius_tube (int or float) – radius of tube (minor radius)
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.TriangularPrism(center, side, height, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
3D Uniform Triangular Prism Axis parallel to z-axis
- Parameters
center (tuple with 3 ints or floats) – center of prism
side (int or float) – side of equilateral base
height (int or float) – height of prism
parameterization (Parameterization) – Parameterization of geometry.
- class modulus.sym.geometry.primitives_3d.VectorizedBoxes(box_bounds, dx=0.0001)[source]
Bases:
Geometry
Vectorized 3D Box/Cuboid for faster surface and interior sampling. This primitive can be used if many boxes are required and is significantly faster then combining many boxes together with Boolean operations.
- Parameters
box_bounds (np.ndarray) – An array specifying the bounds of boxes. Shape of array is [nr_boxes, 3, 2] where the last dim stores the lower and upper bounds respectively.
dx (float) – delta x used for SDF derivative calculations.
Defines base class for all mesh type geometries
- class modulus.sym.geometry.tessellation.Tessellation(mesh, airtight=True, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
Bases:
Geometry
Constructive Tessellation Module that allows sampling on surface and interior of a tessellated geometry.
- Parameters
mesh (Mesh (numpy-stl)) – A mesh that defines the surface of the geometry.
airtight (bool) – If the geometry is airtight or not. If false sample everywhere for interior.
parameterization (Parameterization) – Parameterization of geometry.
- classmethod from_stl(filename, airtight=True, parameterization=<modulus.sym.geometry.parameterization.Parameterization object>)[source]
makes mesh from STL file
- Parameters
filename (str) – filename of mesh.
airtight (bool) – If the geometry is airtight or not. If false sample everywhere for interior.
parameterization (Parameterization) – Parameterization of geometry.