API reference#
- class scifem.PointSource(V: FunctionSpace, points: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]], magnitude: floating | complexfloating = np.float64(1.0))[source]#
Class for defining a point source in a given function space.
- apply_to_vector(b: Function | Vector | Vec, recompute: bool = False)[source]#
Apply the point sources to a vector.
- Parameters:
b – The vector to apply the point sources to.
recompute – If the point sources should be recomputed before applying. Recomputation should only be done if the mesh geometry has been modified.
Note
The user is responsible for forward scattering of the vector after applying the point sources.
Note
If a PETSc vector is passed in, one has to call
b.assemble()
prior to solving the linear system (post scattering).
- scifem.assemble_scalar(J: Form | Form) floating | complexfloating [source]#
Assemble a scalar form and gather result across processes
- Parameters:
form – The form to assemble.
- Returns:
The accumulated value of the assembled form.
- scifem.create_entity_markers(domain: dolfinx.mesh.Mesh, dim: int, entities_list: list[TaggedEntities]) dolfinx.mesh.MeshTags [source]#
Mark entities of specified dimension according to a geometrical marker function.
- Parameters:
domain – A
dolfinx.mesh.Mesh
objectdim – Dimension of the entities to mark
entities_list –
A list of tuples with the following elements:
index 0
: The tag to assign to the entitiesindex 1
: A function that takes a point and returns a boolean array indicating whether the point is inside the entityindex 2
: Optional, if True, the entities will be marked on the boundary
- Returns:
A
dolfinx.mesh.MeshTags
object with the corresponding entities marked. If an entity satisfies multiple input marker functions, it is not deterministic what value the entity gets.
Note:
- scifem.create_real_functionspace(mesh: Mesh, value_shape: tuple[int, ...] = ()) FunctionSpace [source]#
Create a real function space.
- Parameters:
mesh – The mesh the real space is defined on.
value_shape – The shape of the values in the real space.
- Returns:
The real valued function space.
Note
For scalar elements value shape is
()
.
- scifem.dof_to_vertexmap(V: FunctionSpace) ndarray[Any, dtype[int32]] [source]#
Create a map from the degrees of freedom to the vertices of the mesh. As not every degree of freedom is associated with a vertex, every dof that is not associated with a vertex returns -1
- Parameters:
V – The function space
- Returns:
An array mapping local dof i to a local vertex
- scifem.evaluate_function(u: Function, points: Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes], broadcast=True) ndarray[Any, dtype[float64]] [source]#
Evaluate a function at a set of points.
- Parameters:
u – The function to evaluate.
points – The points to evaluate the function at.
broadcast –
If True, the values will be broadcasted to all processes.
Note
Uses a global MPI call to broadcast values, thus this has to be called on all active processes synchronously.
Note
If the function is discontinuous, different processes may return different values for the same point. In this case, the value returned is the maximum value across all processes.
- Returns:
The values of the function evaluated at the points.
- scifem.transfer_meshtags_to_submesh(entity_tag: MeshTags, submesh: Mesh, vertex_to_parent: ndarray[Any, dtype[int32]], cell_to_parent: ndarray[Any, dtype[int32]]) tuple[MeshTags, ndarray[Any, dtype[int32]]] [source]#
Transfer a
entity_tag
from a parent mesh to asubmesh
.- Parameters:
entity_tag – Tag to transfer
submesh – Submesh to transfer tag to
vertex_to_parent – Mapping from submesh vertices to parent mesh vertices
cell_to_parent – Mapping from submesh cells to parent entities
- Returns:
submesh_tag
is the tag on the submesh andsub_to_parent_entity_map
is a mapping from submesh entities in the tag to the corresponding entities in the parent.- Return type:
A tuple (submesh_tag, sub_to_parent_entity_map) where
- scifem.vertex_to_dofmap(V: FunctionSpace) ndarray[Any, dtype[int32]] [source]#
Create a map from the vertices (local to the process) to the correspondning degrees of freedom.
- Parameters:
V – The function space
- Returns:
An array mapping local vertex i to local degree of freedom
Note
If using a blocked space this map is not unrolled for the DofMap block size.