modulus.geometry¶
geometry.geometry¶
Defines base class for all geometries
- class modulus.geometry.geometry.Geometry¶
Bases:
object
Base class for all geometries
geometry.csg.adf¶
- class modulus.geometry.csg.adf.ADF¶
Bases:
torch.nn.modules.module.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[torch.Tensor], radius: float, center: Tuple[float]) torch.Tensor ¶
Computes the pointwise approximate distance for a circle
- points: Tuple[torch.Tensor]
ADF will be computed on these points
- radius: float
Radius of the circle
- center: Tuple[float]
Center of the circle
- omegatorch.Tensor
pointwise approximate distance
- forward(invar)¶
Defines the computation performed at every call.
Should be overridden by all subclasses.
Note
Although the recipe for forward pass needs to be defined within this function, one should call the
Module
instance afterwards instead of this since the former takes care of running the registered hooks while the latter silently ignores them.
- static infinite_line_adf(points: Tuple[torch.Tensor], point_1: Tuple[float], point_2: Tuple[float]) torch.Tensor ¶
Computes the pointwise approximate distance for an infinite line
- 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
- omegatorch.Tensor
pointwise approximate distance
- static line_segment_adf(points: Tuple[torch.Tensor], point_1: Tuple[float], point_2: Tuple[float]) torch.Tensor ¶
Computes the pointwise approximate distance for a line segment
- 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
- omegatorch.Tensor
pointwise approximate distance
- static r_equivalence(omegas: List[torch.Tensor], m: float = 2.0) torch.Tensor ¶
Computes the R-equivalence of a collection of approximate distance functions
- omegasList[torch.Tensor]
List of ADFs used to compute the R-equivalence.
- m: float
Normalization order
- omega_Etorch.Tensor
R-equivalence distance
- static transfinite_interpolation(bases: List[torch.Tensor], indx: int, eps: float = 1e-08) torch.Tensor ¶
Performs transfinite interpolation of the boundary conditions
- 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
- wtorch.Tensor
Interpolation basis corresponding to the input index
- static trimmed_circle_adf(points: Tuple[torch.Tensor], point_1: Tuple[float], point_2: Tuple[float], sign: int, radius: float, center: float) torch.Tensor ¶
Computes the pointwise approximate distance of a trimmed circle
- 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 triming side
- radius: float
Radius of the circle
- center: Tuple[float]
Center of the circle
- omegatorch.Tensor
pointwise approximate distance
geometry.csg.csg¶
Defines base class for all geometries
- class modulus.geometry.csg.csg.ConstructiveSolidGeometry(curves, sdf)¶
Bases:
modulus.geometry.geometry.Geometry
Constructive Solid Geometry Module that allows sampling on surface and interior
- curvelist of Curves
These curves with define the surface or perimeter of geometry.
- sdfSymPy Exprs
SymPy expresion of signed distance function.
- dims()¶
- dimslist of strings
output can be [‘x’], [‘x’,’y’], or [‘x’,’y’,’z’]
- perimeter(param_ranges={}, criteria=None)¶
Perimeter or surface area of geometry.
- param_ranges: dict with of SymPy Symbols and their ranges
If the geometry is parameterized then you can provide ranges for the parameters with this.
- criteriaNone, SymPy boolean exprs
Calculate perimeter discarding areas that don’t satisfy this criteria.
- perimeterfloat
total perimeter or surface area of geometry.
- repeat(spacing, repeat_lower, repeat_higher, center=None)¶
finite repetition of geometry
- spacingfloat, int or SymPy Symbol/Exp
spacing between each repetition.
- repeat_lowertuple of ints
How many repetitions going in negative direction.
- repeat_highertuple of ints
How many repetitions going in positive direction.
- centerNone, list of floats
Do repetition with this center.
- rotate(angle, axis='z', center=None)¶
rotates geometry
- anglefloat, int, SymPy Symbol/Exp
Angle to rotate geometry.
- axisstr
Rotate around this axis.
- centerNone, list of floats
Do rotation around this center
- sample_boundary(nr_points, criteria=None, param_ranges={}, quasirandom=False)¶
Samples the surface or perimeter of the geometry.
- nr_pointsint
number of points to sample on boundary.
- criteriaNone, SymPy boolean exprs
Only sample points that satisfy this criteria.
- param_ranges: dict with of SymPy Symbols and their ranges
If the geometry is parameterized then you can provide ranges for the parameters with this.
- quasirandombool
If true then sample the points using the Halton sequences. Default is False.
- pointsdictionary of numpy arrays
Dictionary contain a point cloud sample 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’])
- sample_interior(nr_points, bounds, criteria=None, param_ranges={}, quasirandom=False)¶
Samples the interior of the geometry.
- nr_pointsint
number of points to sample.
- boundsdict
bounds to sample points from. For example, bounds = {‘x’: (0, 1), ‘y’: (0, 1)}.
- criteriaNone, SymPy boolean exprs
Only sample points that satisfy this criteria.
- param_ranges: dict with of SymPy Symbols and their ranges
If the geometry is parameterized then you can provide ranges for the parameters with this.
- quasirandombool
If true then sample the points using the Halton sequences. Default is False.
- pointsdictionary of numpy arrays
Dictionary contain a point cloud sample 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’])
- scale(x, center=None)¶
scale geometry
- xfloat, int, SymPy Symbol/Exp
Scale factor.
- centerNone, list of floats
Do scaling with around this center.
- translate(xyz)¶
translate geometry
- xyztuple of float, int or SymPy Symbol/Exp
translate geometry by these values.
geometry.csg.csg_1d¶
- class modulus.geometry.csg.csg_1d.Line1D(point_1, point_2)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
1D Line along x-axis
- point_1int or float
lower bound point of line
- point_2int or float
upper bound point of line
- class modulus.geometry.csg.csg_1d.Point1D(point)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
1D Point along x-axis
- pointint or float
x coordinate of the point
geometry.csg.csg_2d¶
Primitives for 2D geometries see https://www.iquilezles.org/www/articles/distfunctions/distfunctions.html
- class modulus.geometry.csg.csg_2d.Channel2D(point_1, point_2)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Channel (no bounding curves in x-direction)
- point_1tuple with 2 ints or floats
lower bound point of channel
- point_2tuple with 2 ints or floats
upper bound point of channel
- class modulus.geometry.csg.csg_2d.Circle(center, radius)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Circle
- centertuple with 2 ints or floats
center point of circle
- radiusint or float
radius of circle
- class modulus.geometry.csg.csg_2d.Ellipse(center, major, minor)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Ellipse
- centertuple with 2 ints or floats
center point of circle
- radiusint or float
radius of circle
- class modulus.geometry.csg.csg_2d.Line(point_1, point_2, normal)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Line parallel to y-axis
- point_1tuple with 2 ints or floats
lower bound point of line segment
- point_2tuple with 2 ints or floats
upper bound point of line segment
- normalint or float
normal direction of line (+1 or -1)
- class modulus.geometry.csg.csg_2d.Rectangle(point_1, point_2)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Rectangle
- point_1tuple with 2 ints or floats
lower bound point of rectangle
- point_2tuple with 2 ints or floats
upper bound point of rectangle
- class modulus.geometry.csg.csg_2d.Triangle(center, base, height)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Isosceles Triangle Symmetrical axis parallel to y-axis
- centertuple with 2 ints or floats
center of base of triangle
- baseint or float
base of triangle
- heightint or float
height of triangle
geometry.csg.csg_3d¶
Primitives for 3D geometries see https://www.iquilezles.org/www/articles/distfunctions/distfunctions.html
- class modulus.geometry.csg.csg_3d.Box(point_1, point_2)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Box/Cuboid
- point_1tuple with 3 ints or floats
lower bound point of box
- point_2tuple with 3 ints or floats
upper bound point of box
- class modulus.geometry.csg.csg_3d.Channel(point_1, point_2)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Channel (no bounding surfaces in x-direction)
- point_1tuple with 3 ints or floats
lower bound point of channel
- point_2tuple with 3 ints or floats
upper bound point of channel
- class modulus.geometry.csg.csg_3d.Cone(center, radius, height)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Cone Axis parallel to z-axis
- centertuple with 3 ints or floats
base center of cone
- radiusint or float
base radius of cone
- heightint or float
height of cone
- class modulus.geometry.csg.csg_3d.Cylinder(center, radius, height)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Cylinder Axis parallel to z-axis
- centertuple with 3 ints or floats
center of cylinder
- radiusint or float
radius of cylinder
- heightint or float
height of cylinder
- class modulus.geometry.csg.csg_3d.ElliCylinder(center, a, b, height)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
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
- centertuple with 3 ints or floats
center of base of ellipse
- aint or float
semi-major axis of ellipse
- bint or float
semi-minor axis of ellipse
- heightint or float
height of elliptical cylinder
- class modulus.geometry.csg.csg_3d.IsoTriangularPrism(center, base, height, height_prism)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
2D Isosceles Triangular Prism Symmetrical axis parallel to y-axis
- centertuple with 3 ints or floats
center of base of triangle
- baseint or float
base of triangle
- heightint or float
height of triangle
- height_prismint or float
height of triangular prism
- class modulus.geometry.csg.csg_3d.Plane(point_1, point_2, normal)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Plane perpendicular to x-axis
- point_1tuple with 3 ints or floats
lower bound point of plane
- point_2tuple with 3 ints or floats
upper bound point of plane
- normalint or float
normal direction of plane (+1 or -1)
- class modulus.geometry.csg.csg_3d.Sphere(center, radius)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Sphere
- centertuple with 3 ints or floats
center of sphere
- radiusint or float
radius of sphere
- class modulus.geometry.csg.csg_3d.Tetrahedron(center, radius)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
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
- centertuple with 3 ints or floats
centroid of tetrahedron
- radiusint or float
radius of circumscribed sphere
- class modulus.geometry.csg.csg_3d.Torus(center, radius, radius_tube)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Torus
- centertuple with 3 ints or floats
center of torus
- radiusint or float
distance from center to center of tube (major radius)
- radius_tubeint or float
radius of tube (minor radius)
- class modulus.geometry.csg.csg_3d.TriangularPrism(center, side, height)¶
Bases:
modulus.geometry.csg.csg.ConstructiveSolidGeometry
3D Uniform Triangular Prism Axis parallel to z-axis
- centertuple with 3 ints or floats
center of prism
- sideint or float
side of equilateral base
- heightint or float
height of prism
geometry.csg.curves¶
Defines parameterized curves in SymPy
- class modulus.geometry.csg.curves.Curve(functions, ranges, area, criteria=True)¶
Bases:
object
Parameterized curves that keeps track of normals and area/perimeter
- functionsdictionary of SymPy Exprs
Parameterized curve in 1, 2 or 3 dimensions. For example, a circle might have:
functions = {'x': cos(theta), 'y': sin(theta), 'normal_x': cos(theta), 'normal_y': sin(theta)}
TODO: refactor to remove normals.
- rangesdictionary 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)}.
- areafloat, int, SymPy Exprs
The surface area/perimeter of the curve.
- criteriaSymPy 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.
- approx_area(param_ranges={}, criteria=None, approx_nr=10000)¶
- param_ranges: dict with of SymPy Symbols and their ranges
If the curve is parameterized then you can provide ranges for the parameters with this.
- criteriaNone, SymPy boolean exprs
Calculate area discarding regions that don’t satisfy this criteria.
- approx_nrint
Area might be difficult to compute if parameterized. In this case the area is approximated by sampleing self.area, approx_nr number of times. This amounts to monte carlo integration.
- areafloat
area of curve
- dims()¶
- dimslist of strings
output can be [‘x’], [‘x’,’y’], or [‘x’,’y’,’z’]
- rotate(angle, axis='z')¶
rotate curve
- anglefloat, SymPy Symbol/Exprs
angle of rotation
- axisstr
axis of rotation
- sample(nr_points, param_ranges={}, quasirandom=False)¶
sample curve
- nr_pointsint
Number of points per area to sample.
- param_ranges: dict with of SymPy Symbols and their ranges
If the curve is parameterized then you can give the ranges for these parameters with this.
- quasirandombool
If true then sample the points using the Halton sequences. Default is False.
- scale(x)¶
scale curve
- xfloat, SymPy Symbol/Exprs
scale factor.
- translate(xyz)¶
translate curve
- xyztuple of floats, ints, SymPy Symbol/Exprs
translate curve by these values.
geometry.tessellation.tessellation¶
Defines base class for all mesh type geometries
- class modulus.geometry.tessellation.tessellation.Tessellation(mesh, airtight=True, epsilon=1e-06)¶
Bases:
modulus.geometry.geometry.Geometry
Constructive Tessellation Module that allows sampling on surface and interior of a tessellated geometry.
- meshMesh (numpy-stl)
A mesh that defines the surface of the geometry.
- airtightbool
If the geometry is airtight or not. If false sample everywere for interior.
- epsilonfloat
Minimum distance to wall when sampling interior of mesh
- bounds()¶
x, y, z bounds of geometry.
- flip(axis)¶
flips geometry
- axisstr
‘x’, ‘y’, or ‘z’ for which axis to flip. If its ‘x’ then the geometry will be flipped so x_new = -x
- classmethod from_stl(filename, airtight=True, epsilon=1e-06)¶
makes mesh from STL file
- filenamestr
filename of mesh.
- airtightbool
If the geometry is airtight or not. If false sample everywere for interior.
- epsilonfloat
Minimum distance to wall when sampling interior of mesh
- rotate(vector, radians)¶
rotates geometry
- vectortuple of float, int
vector of rotation.
- radiansfloat, int
rotation in radians.
- sample_boundary(nr_points, criteria=None, param_ranges={}, quasirandom=False)¶
Samples the surface of the geometry.
- nr_pointsint
number of points to sample.
- criteriaNone, SymPy Boolean expressions
Only sample points that satisfy this criteria.
- param_ranges: dict with of SymPy Symbols and their ranges
If the geometry is parameterized then you can provide ranges for the parameters with this.
- quasirandom: bool
if true use quasirandom sampling (TODO not implemented)
- pointsdictionary of numpy arrays
Dictionary contain a point cloud sample uniformly. For example:
points = {'x': np.ndarray (N, 1), 'y': np.ndarray (N, 1), 'z': np.ndarray (N, 1), 'normal_x': np.ndarray (N, 1), 'normal_y': np.ndarray (N, 1), 'normal_z': 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’])
- sample_interior(nr_points, bounds=None, criteria=None, compute_distance_field=True, param_ranges={}, quasirandom=False)¶
Samples the interior of the geometry.
- nr_pointsint
number of points to sample.
- boundsdict
bounds to sample points from. For example, bounds = {‘x’: (0, 1), ‘y’: (0, 1)}.
- criteriaNone, SymPy boolean exprs
Only sample points that satisfy this criteria.
- compute_distance_fieldbool
If true then return the signed distance field and its derivatives along with the point cloud.
- param_ranges: dict with of SymPy Symbols and their ranges
If the geometry is parameterized then you can provide ranges for the parameters with this.
- quasirandombool
If true then sample the points using the Halton sequences. Default is False.
- pointsdictionary of numpy arrays
Dictionary contain a point cloud sample uniformly. For example:
points = {'x': np.ndarray (N, 1), 'y': np.ndarray (N, 1), 'z': np.ndarray (N, 1), 'sdf': np.ndarray (N, 1), 'sdf__x': np.ndarray (N, 1), # derivatives of SDF 'sdf__y': np.ndarray (N, 1), 'sdf__z': np.ndarray (N, 1), 'normal_x': np.ndarray (N, 1), 'normal_y': np.ndarray (N, 1), 'normal_z': 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’])
- scale(factor, center=None)¶
scale geometry
- factorfloat, int, SymPy Symbol/Exp
Scale factor.
- centerNone, list of floats
Do scaling with around this center.
- surface_area()¶
Perimeter or surface area of geometry.
- surface_areafloat
total surface area of geometry.
- translate(xyz)¶
translates geometry
- xyztuple of float, int
translate geometry by these values.