symfem.references
¶
Reference elements.
Module Contents¶
Classes¶
A reference cell on which a finite element can be defined. |
|
A point. |
|
An interval. |
|
A triangle. |
|
A tetrahedron. |
|
A quadrilateral. |
|
A hexahedron. |
|
A (triangular) prism. |
|
A (square-based) pyramid. |
|
A polygon on a barycentric dual grid. |
Functions¶
|
Check which side of a line or plane a set of points are. |
|
Subtract. |
|
Add. |
|
Compute the dot product. |
|
Compute the cross product. |
|
Find the norm of a vector. |
|
Normalise a vector. |
Attributes¶
- symfem.references.LatticeWithLines¶
- symfem.references.IntLimits¶
- exception symfem.references.NonDefaultReferenceError¶
Bases:
NotImplementedError
Exception to be thrown when an element can only be created on the default reference.
- symfem.references._which_side(vs: symfem.geometry.SetOfPoints, p: symfem.geometry.PointType, q: symfem.geometry.PointType) int | None ¶
Check which side of a line or plane a set of points are.
- Parameters:
vs – The set of points
p – A point on the line or plane
q – Another point on the line (2D) or the normal to the plane (3D)
- Returns:
2 if the points are all to the left, 1 if the points are all to the left or on the line, 0 if the points are all on the line, -1 if the points are all to the right or on the line, -1 if the points are all to the right, None if there are some points on either side.
- symfem.references._vsub(v: symfem.geometry.PointTypeInput, w: symfem.geometry.PointTypeInput) symfem.geometry.PointType ¶
Subtract.
- Parameters:
v – A vector
w – A vector
- Returns:
The vector v - w
- symfem.references._vadd(v: symfem.geometry.PointTypeInput, w: symfem.geometry.PointTypeInput) symfem.geometry.PointType ¶
Add.
- Parameters:
v – A vector
w – A vector
- Returns:
The vector v + w
- symfem.references._vdot(v: symfem.geometry.PointTypeInput, w: symfem.geometry.PointTypeInput) sympy.core.expr.Expr ¶
Compute the dot product.
- Parameters:
v – A vector
w – A vector
- Returns:
The scalar v.w
- symfem.references._vcross(v_in: symfem.geometry.PointTypeInput, w_in: symfem.geometry.PointTypeInput) symfem.geometry.PointType ¶
Compute the cross product.
- Parameters:
v – A vector
w – A vector
- Returns:
The vector v x w
- symfem.references._vnorm(v_in: symfem.geometry.PointTypeInput) sympy.core.expr.Expr ¶
Find the norm of a vector.
- Parameters:
v_in – A vector
- Returns:
The norm of v_in
- symfem.references._vnormalise(v_in: symfem.geometry.PointTypeInput) symfem.geometry.PointType ¶
Normalise a vector.
- Parameters:
v_in – A vector
- Returns:
A unit vector pointing in the same direction as v_in
- class symfem.references.Reference(vertices: symfem.geometry.SetOfPointsInput = ())¶
Bases:
abc.ABC
A reference cell on which a finite element can be defined.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- __build__(tdim: int, name: str, origin: symfem.geometry.PointTypeInput, axes: symfem.geometry.SetOfPointsInput, reference_vertices: symfem.geometry.SetOfPointsInput, vertices: symfem.geometry.SetOfPointsInput, edges: Tuple[Tuple[int, int], Ellipsis], faces: Tuple[Tuple[int, Ellipsis], Ellipsis], volumes: Tuple[Tuple[int, Ellipsis], Ellipsis], sub_entity_types: List[List[str] | str | None], simplex: bool = False, tp: bool = False)¶
Create a reference cell.
- Parameters:
tdim – The topological dimension of the cell
name – The name of the cell
origin – The coordinates of the origin
axes – Vectors representing the axes of the cell
reference_vertices – The vertices of the default version of this cell
vertices – The vertices of this cell
edges – Pairs of vertex numbers that form the edges of the cell
faces – Tuples of vertex numbers that form the faces of the cell
volumes – Tuples of vertex numbers that form the volumes of the cell
sub_entity_types – The cell types of each sub-entity of the cell
simplex – Is the cell a simplex (interval/triangle/tetrahedron)?
tp – Is the cell a tensor product (interval/quadrilateral/hexahedron)?
- __eq__(other: object) bool ¶
Check if two references are equal.
- __hash__() int ¶
Check if two references are equal.
- intersection(other: Reference) Reference | None ¶
Get the intersection of two references.
- Returns:
A reference element that is the intersection
- abstract default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- abstract make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- make_lattice_float(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- make_lattice_with_lines_float(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- z_ordered_entities() List[List[Tuple[int, int]]] ¶
Get the subentities of the cell in back-to-front plotting order.
- Returns:
List of lists of subentity dimensions and numbers
- z_ordered_entities_extra_dim() List[List[Tuple[int, int]]] ¶
Get the subentities in back-to-front plotting order when using an extra dimension.
- Returns:
List of lists of subentity dimensions and numbers
- get_point(reference_coords: symfem.geometry.PointType) Tuple[sympy.core.expr.Expr, Ellipsis] ¶
Get a point in the reference from reference coordinates.
- Parameters:
reference_coords – The reference coordinates
- Returns:
A point in the cell
- abstract integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- abstract get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- abstract get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- get_map_to_self() symfem.geometry.PointType ¶
Get the map from the canonical reference to this reference.
- Returns:
The map
- abstract _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- get_inverse_map_to_self() symfem.geometry.PointType ¶
Get the inverse map from the canonical reference to this reference.
- Returns:
The map
- abstract _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- abstract volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- midpoint() symfem.geometry.PointType ¶
Calculate the midpoint.
- Returns:
The midpoint of the cell
- jacobian() sympy.core.expr.Expr ¶
Calculate the Jacobian.
- Returns:
The Jacobian
- scaled_axes() symfem.geometry.SetOfPoints ¶
Return the unit axes of the reference.
- Returns:
The axes
- tangent() symfem.geometry.PointType ¶
Calculate the tangent to the element.
- Returns:
The tangent
- normal() symfem.geometry.PointType ¶
Calculate the normal to the element.
- Returns:
The normal
- sub_entities(dim: int | None = None, codim: int | None = None) Tuple[Tuple[int, Ellipsis], Ellipsis] ¶
Get the sub-entities of a given dimension.
- Parameters:
dim – The dimension of the sub-entity
codim – The co-dimension of the sub-entity
- Returns:
A tuple of tuples of vertex numbers
- sub_entity_count(dim: int | None = None, codim: int | None = None) int ¶
Get the number of sub-entities of a given dimension.
- Parameters:
dim – the dimension of the sub-entity
codim – the codimension of the sub-entity
- Returns:
The number of sub-entities
- sub_entity(dim: int, n: int, reference_vertices: bool = False) Any ¶
Get the sub-entity of a given dimension and number.
- Parameters:
dim – the dimension of the sub-entity
n – The sub-entity number
reference_vertices – Should the reference vertices be used?
- Returns:
The sub-entity
- at_vertex(point: symfem.geometry.PointType) bool ¶
Check if a point is a vertex of the reference.
- Parameters:
point – The point
- Returns:
Is the point a vertex?
- on_edge(point_in: symfem.geometry.PointType) bool ¶
Check if a point is on an edge of the reference.
- Parameters:
point_in – The point
- Returns:
Is the point on an edge?
- on_face(point_in: symfem.geometry.PointType) bool ¶
Check if a point is on a face of the reference.
- Parameters:
point_in – The point
- Returns:
Is the point on a face?
- abstract contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- map_polyset_from_default(poly: List[symfem.functions.FunctionInput]) List[symfem.functions.FunctionInput] ¶
Map the polynomials from the default reference element to this reference.
- plot_entity_diagrams(filename: str | List[str], plot_options: Dict[str, Any] = {}, **kwargs: Any)¶
Plot diagrams showing the entity numbering of the reference.
- class symfem.references.Point(vertices: symfem.geometry.SetOfPointsInput = ((),))¶
Bases:
Reference
A point.
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- abstract make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Interval(vertices: symfem.geometry.SetOfPointsInput = ((0,), (1,)))¶
Bases:
Reference
An interval.
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Triangle(vertices: symfem.geometry.SetOfPointsInput = ((0, 0), (1, 0), (0, 1)))¶
Bases:
Reference
A triangle.
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- z_ordered_entities_extra_dim() List[List[Tuple[int, int]]] ¶
Get the subentities in back-to-front plotting order when using an extra dimension.
- Returns:
List of lists of subentity dimensions and numbers
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Tetrahedron(vertices: symfem.geometry.SetOfPointsInput = ((0, 0, 0), (1, 0, 0), (0, 1, 0), (0, 0, 1)))¶
Bases:
Reference
A tetrahedron.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- z_ordered_entities() List[List[Tuple[int, int]]] ¶
Get the subentities of the cell in back-to-front plotting order.
- Returns:
List of lists of subentity dimensions and numbers
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Quadrilateral(vertices: symfem.geometry.SetOfPointsInput = ((0, 0), (1, 0), (0, 1), (1, 1)))¶
Bases:
Reference
A quadrilateral.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- z_ordered_entities_extra_dim() List[List[Tuple[int, int]]] ¶
Get the subentities in back-to-front plotting order when using an extra dimension.
- Returns:
List of lists of subentity dimensions and numbers
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Hexahedron(vertices: symfem.geometry.SetOfPointsInput = ((0, 0, 0), (1, 0, 0), (0, 1, 0), (1, 1, 0), (0, 0, 1), (1, 0, 1), (0, 1, 1), (1, 1, 1)))¶
Bases:
Reference
A hexahedron.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- z_ordered_entities() List[List[Tuple[int, int]]] ¶
Get the subentities of the cell in back-to-front plotting order.
- Returns:
List of lists of subentity dimensions and numbers
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Prism(vertices: symfem.geometry.SetOfPointsInput = ((0, 0, 0), (1, 0, 0), (0, 1, 0), (0, 0, 1), (1, 0, 1), (0, 1, 1)))¶
Bases:
Reference
A (triangular) prism.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- z_ordered_entities() List[List[Tuple[int, int]]] ¶
Get the subentities of the cell in back-to-front plotting order.
- Returns:
List of lists of subentity dimensions and numbers
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.Pyramid(vertices: symfem.geometry.SetOfPointsInput = ((0, 0, 0), (1, 0, 0), (0, 1, 0), (1, 1, 0), (0, 0, 1)))¶
Bases:
Reference
A (square-based) pyramid.
- property clockwise_vertices: symfem.geometry.SetOfPoints¶
Get list of vertices in clockwise order.
- Returns:
A list of vertices
- z_ordered_entities() List[List[Tuple[int, int]]] ¶
Get the subentities of the cell in back-to-front plotting order.
- Returns:
List of lists of subentity dimensions and numbers
- default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- abstract make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- class symfem.references.DualPolygon(number_of_triangles: int, vertices: symfem.geometry.SetOfPointsInput | None = None)¶
Bases:
Reference
A polygon on a barycentric dual grid.
- reference_origin: symfem.geometry.PointType¶
- abstract contains(point: symfem.geometry.PointType) bool ¶
Check if a point is contained in the reference.
- Parameters:
point – The point
- Returns:
Is the point contained in the reference?
- abstract integration_limits(vars: symfem.symbols.AxisVariablesNotSingle = t) IntLimits ¶
Get the limits for an integral over this reference.
- Parameters:
vars – The variables to use for each direction
- Returns:
Integration limits that can be passed into sympy.integrate
- abstract get_map_to(vertices: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the map from the reference to a cell.
- Parameters:
vertices – The vertices of a call
- Returns:
The map
- abstract get_inverse_map_to(vertices_in: symfem.geometry.SetOfPointsInput) symfem.geometry.PointType ¶
Get the inverse map from a cell to the reference.
- Parameters:
vertices_in – The vertices of a cell
- Returns:
The inverse map
- abstract _compute_map_to_self() symfem.geometry.PointType ¶
Compute the map from the canonical reference to this reference.
- Returns:
The map
- abstract _compute_inverse_map_to_self() symfem.geometry.PointType ¶
Compute the inverse map from the canonical reference to this reference.
- Returns:
The map
- abstract volume() sympy.core.expr.Expr ¶
Calculate the volume.
- Returns:
The volume of the cell
- abstract default_reference() Reference ¶
Get the default reference for this cell type.
- Returns:
The default reference
- make_lattice(n: int) symfem.geometry.SetOfPoints ¶
Make a lattice of points.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points offset from the edge of the cell
- abstract make_lattice_with_lines(n: int) LatticeWithLines ¶
Make a lattice of points, and a list of lines connecting them.
- Parameters:
n – The number of points along each edge
- Returns:
A lattice of points including the edges of the cell Pairs of point numbers that make a mesh of lines across the cell
- z_ordered_entities_extra_dim() List[List[Tuple[int, int]]] ¶
Get the subentities in back-to-front plotting order when using an extra dimension.
- Returns:
List of lists of subentity dimensions and numbers