API Reference#
- class oasisx.DirichletBC(value: float64 | Constant | Callable[[ndarray[Any, dtype[float64]]], ndarray[Any, dtype[float64]]], method: LocatorMethod, marker: Tuple[MeshTags, int32] | Callable[[ndarray[Any, dtype[float64]]], ndarray[Any, dtype[bool]]])[source]#
Create a Dirichlet boundary condition based on topological or geometrical info from the mesh This boundary condtion should only be used for velocity function spaces.
- Parameters:
value – The value the degrees of freedom should have. It can be a float, a dolfinx.fem.Constant or a lambda-function.
method – Locator method for marker.
marker – If
oasisx.LocatorMethod.TOPOLOGICAL
the input should be a tuple of a mesh tag and the corresponding value for the entities to assign Dirichlet conditions to. Ifoasisx.LocatorMethod.GEOMETRICAL
the input is a lambda function taking in x, y, and z coordinates (ordered as [[x0, x1,…,xn], [y0,…,yn], [z0,…,zn]]), and return a boolean marker where the ith entry indicates if [xi, yi, zi] should have a diriclet bc.
Example
Assigning a topological condition
entities = np.array([0,3,8],dtype=np.int32) values = np.full_like(entities, 2, dtype=np.int32) mt = dolfinx.fem.meshtags(mesh, mesh.topology.dim-1, entities, values) bc = DirichletBC(5., LocatorMethod.TOPOLOGICAL, (mt, 2))
Assigning a geometrical condition
bc = DirchletBC(dolfinx.fem.Constant(mesh, 3.), lambda x: np.isclose(x[0]))
- class oasisx.FractionalStep_AB_CN(mesh: Mesh, u_element: Tuple[str, int] | FiniteElement, p_element: Tuple[str, int] | FiniteElement, bcs_u: List[List[DirichletBC]], bcs_p: List[PressureBC], rotational: bool = False, solver_options: dict | None = None, jit_options: dict | None = None, body_force: Expr | None = None, options: dict | None = None)[source]#
Create the fractional step solver with Adam-Bashforth linearization of the convective term, and Crank-Nicholson time discretization.
- Parameters:
mesh – The computational domain
u_element – A tuple describing the finite element family and degree of the element used for the velocity
p_element – A tuple describing the finite element family and degree of the element used for the pressure
bcs_u – List of Dirichlet BCs for each component of the velocity
bcs_p – List of Dirichlet BCs for the pressure
rotational – If True, use rotational form of pressure update
solver_options – Dictionary with keys ‘tentative’, ‘pressure’ and ‘scalar’, where each key leads to a dictionary of of PETSc options for each problem
jit_options – Just in time parameters, to optimize form compilation
options – Options for the Oasis solver. “low_memory_version” True/False changes if \(\int\nabla_k p^* v~\mathrm{d}x\) is assembled as True: directly into a vector or False: matrix-vector product. Default value: True ‘“bc_topological”` True/False. changes how the Dirichlet dofs are located. If True facet_markers has to be supplied.
body_force – List of forces acting in each direction (x,y,z)
- assemble_first(dt: float, nu: float)[source]#
Assembly routine for first iteration of pressure/velocity system.
\[A=\frac{1}{\delta t}M + \frac{1}{2} C +\frac{1}{2}\nu K\]where M is the mass matrix, K the stiffness matrix and C the convective term. We also assemble parts of the right hand side:
\[b_k=\left(\frac{1}{\delta t}M 0 \frac{1}{2} C -\frac{1}{2}\nu K\right)u_k^{n-1} + \int_{\Omega} f_k^{n-\frac{1}{2}}v_k~\mathrm{d}x + \int h^{n-\frac{1}{2}} n_k \nabla_k v~\mathrm{d}x\]- Parameters:
dt – The time step
nu – The kinematic viscosity
- pressure_assemble(dt: float)[source]#
Assemble RHS for pressure correction term
\[b_c = \int_{\Omega} \mathrm{div} u^l q~\mathrm{d}x.\]
- pressure_solve(nu: float | None = None, rotational: bool = False) int32 [source]#
Solve pressure correction problem
- solve(dt: float, nu: float, max_error: float = 1e-12, max_iter: int = 10)[source]#
Propagate splitting scheme one time step
- Parameters:
dt – The time step
nu – The kinematic velocity
max_error – The maximal difference for inner iterations of solving u
max_iter – Maximal number of inner iterations for u
- property u#
Return the solution to the tentative velocity equation as a vector function
- velocity_tentative_assemble()[source]#
Assemble RHS of tentative velocity equation computing the \(p^*\)-dependent term of
\[b_k \mathrel{+}= b(u_k^{n-1}) + \int_{\Omega} p^*\frac{\partial v_k}{\partial x_k}~\mathrm{d}x.\]\(b(u_k^{n-1})\) is computed in
FractionalStep_AB_CN.assemble_first
.
- class oasisx.LocatorMethod(value, names=<not given>, *values, module=None, qualname=None, type=None, start=1, boundary=None)[source]#
Search methods for Dirichlet BCs
- GEOMETRICAL = 1#
- TOPOLOGICAL = 2#
- class oasisx.PressureBC(value: float64 | Constant | Callable[[ndarray[Any, dtype[float64]]], ndarray[Any, dtype[float64]]], marker: Tuple[MeshTags, int32 | Tuple[int32]])[source]#
Create a Pressure boundary condition (natural bc) based on a set of facets of a mesh and some value.
- Parameters:
value – The value the degrees of freedom should have. It can be a float, a dolfinx.fem.Constant or a lambda-function. If value is a lambda-function it is interpolated into the pressure space.
marker – Tuple of a mesh tag and the corresponding value for the entities to assign Dirichlet conditions to. The meshtag dimension has to be mesh.topology.dim -1.
Example
Assigning a constant condition
entities = np.array([0,3,8],dtype=np.int32) values = np.full_like(entities, 2, dtype=np.int32) mt = dolfinx.fem.meshtags(mesh, mesh.topology.dim-1, entities, values) bc = DirichletBC(5., (mt, 2))
Assigning a time-dependent condition
class OutletPressure(): def __init__(self, t:float): self.t = t def eval(self, x: numpy.typing.NDArray[np.float64]): return self.t*x[0] p = OutletPressure(0.) entities = np.array([0,3,8],dtype=np.int32) values = np.full_like(entities, 2, dtype=np.int32) mt = dolfinx.fem.meshtags(mesh, mesh.topology.dim-1, entities, values) bc = PressureBC(p, (mt, 2)) Q = dolfinx.fem.functionspace(mesh, ("Lagrange", 1)) # Create appropriate structures for assigning bcs bc.create_bc(Q) # Update time in bc p.t = 1
- class oasisx.Projector(function: Expr, space: FunctionSpace, bcs: List[DirichletBC] | None = None, petsc_options: dict | None = None, jit_options: dict | None = None, form_compiler_options: dict | None = None, metadata: dict | None = None)[source]#
Projector for a given function. Solves Ax=b, where
u, v = ufl.TrialFunction(Space), ufl.TestFunction(space) dx = ufl.Measure("dx", metadata=metadata) A = inner(u, v) * dx b = inner(function, v) * dx
- Parameters:
function – UFL expression of function to project
space – Space to project function into
bcs – List of BCS to apply to projection
petsc_options – Options to pass to PETSc
jit_options – Options to pass to just in time compiler
form_compiler_options – Options to pass to the form compiler
metadata – Data to pass to the integration measure