dolfinx_mpc#

Main module for DOLFINX_MPC

class dolfinx_mpc.LinearProblem(a: Form, L: Form, mpc: MultiPointConstraint, bcs: List[DirichletBC] | None = None, u: Function | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]#

Class for solving a linear variational problem with multi point constraints of the form a(u, v) = L(v) for all v using PETSc as a linear algebra backend.

Parameters:
  • a – A bilinear UFL form, the left hand side of the variational problem.

  • L – A linear UFL form, the right hand side of the variational problem.

  • mpc – The multi point constraint.

  • bcs – A list of Dirichlet boundary conditions.

  • u

    The solution function. It will be created if not provided. The function has to be based on the functionspace in the mpc, i.e.

    u = dolfinx.fem.Function(mpc.function_space)
    

  • petsc_options – Parameters that is passed to the linear algebra backend PETSc. #type: ignore For available choices for the ‘petsc_options’ kwarg, see the PETSc-documentation https://www.mcs.anl.gov/petsc/documentation/index.html.

  • form_compiler_options – Parameters used in FFCx compilation of this form. Run ffcx –help at the commandline to see all available options. Takes priority over all other parameter values, except for scalar_type which is determined by DOLFINx.

  • jit_options – Parameters used in CFFI JIT compilation of C code generated by FFCx. See FEniCS/dolfinx for all available parameters. Takes priority over all other parameter values.

Examples

Example usage:

problem = LinearProblem(a, L, mpc, [bc0, bc1],
                        petsc_options={"ksp_type": "preonly", "pc_type": "lu"})
solve() Function[source]#

Solve the problem.

Returns:

Function containing the solution

class dolfinx_mpc.MultiPointConstraint(V: ~dolfinx.fem.function.FunctionSpace, dtype: ~numpy.floating = <class 'numpy.float64'>)[source]#

Hold data for multi point constraint relation ships, including new index maps for local assembly of matrices and vectors.

Parameters:
  • V – The function space

  • dtype – The dtype of the underlying functions

add_constraint(V: FunctionSpace, slaves: ndarray[Any, dtype[int32]], masters: ndarray[Any, dtype[int64]], coeffs: ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]] | ndarray[Any, dtype[complex64]] | ndarray[Any, dtype[complex128]], owners: ndarray[Any, dtype[int32]], offsets: ndarray[Any, dtype[int32]])[source]#

Add new constraint given by numpy arrays.

Parameters:
  • V – The function space for the constraint

  • slaves – List of all slave dofs (using local dof numbering) on this process

  • masters – List of all master dofs (using global dof numbering) on this process

  • coeffs – The coefficients corresponding to each master.

  • owners – The process each master is owned by.

  • offsets

    Array indicating the location in the masters array for the i-th slave in the slaves arrays, i.e.

    masters_of_owned_slave[i] = masters[offsets[i]:offsets[i+1]]
    

add_constraint_from_mpc_data(V: FunctionSpace, mpc_data: mpc_data_double | mpc_data_float | mpc_data_complex_double | mpc_data_complex_float | MPCData)[source]#

Add new constraint given by an dolfinc_mpc.cpp.mpc.mpc_data-object

backsubstitution(u: Function | Vec) None[source]#

For a Function, impose the multi-point constraint by backsubstiution. This function is used after solving the reduced problem to obtain the values at the slave degrees of freedom

Note

It is the users responsibility to destroy the PETSc vector

Parameters:

u – The input function

property cell_to_slaves#

Returns an dolfinx.cpp.graph.AdjacencyList_int32 whose ith node corresponds to the ith cell (local to process), and links the corresponding slave degrees of freedom in the cell (local to process).

Examples

cell_to_slaves = mpc.cell_to_slaves()
slaves_in_cell_i = cell_to_slaves.links(i)
coefficients() ndarray[Any, dtype[float32]] | ndarray[Any, dtype[float64]] | ndarray[Any, dtype[complex64]] | ndarray[Any, dtype[complex128]][source]#

Returns a vector containing the coefficients for the constraint, and the corresponding offsets for the ith degree of freedom.

Examples

coeffs, offsets = mpc.coefficients()
coeffs_of_slave_i = coeffs[offsets[i]:offsets[i+1]]
create_contact_inelastic_condition(meshtags: MeshTags_int32, slave_marker: int, master_marker: int, eps2: float = 1e-20)[source]#

Create a contact inelastic condition between two sets of facets marker with individual markers. The interfaces should be within machine precision of eachother, but the vertices does not need to align. The condition created is \(u_s = u_m\) where s is the restriction to the slave facets, m to the master facets.

Parameters:
  • meshtags – The meshtags of the set of facets to tie together

  • slave_marker – The marker of the slave facets

  • master_marker – The marker of the master facets

  • eps2 – The tolerance for the squared distance between cells to be considered as a collision

create_contact_slip_condition(meshtags: MeshTags, slave_marker: int, master_marker: int, normal: Function, eps2: float = 1e-20)[source]#

Create a slip condition between two sets of facets marker with individual markers. The interfaces should be within machine precision of eachother, but the vertices does not need to align. The condition created is \(u_s \cdot normal_s = u_m \cdot normal_m\) where s is the restriction to the slave facets, m to the master facets.

Parameters:
  • meshtags – The meshtags of the set of facets to tie together

  • slave_marker – The marker of the slave facets

  • master_marker – The marker of the master facets

  • normal – The function used in the dot-product of the constraint

  • eps2 – The tolerance for the squared distance between cells to be considered as a collision

create_general_constraint(slave_master_dict: Dict[bytes, Dict[bytes, float]], subspace_slave: int | None = None, subspace_master: int | None = None)[source]#
Parameters:
  • V – The function space

  • slave_master_dict – Nested dictionary, where the first key is the bit representing the slave dof’s coordinate in the mesh. The item of this key is a dictionary, where each key of this dictionary is the bit representation of the master dof’s coordinate, and the item the coefficient for the MPC equation.

  • subspace_slave – If using mixed or vector space, and only want to use dofs from a sub space as slave add index here

  • subspace_master – Subspace index for mixed or vector spaces

Example

If the dof D located at [d0, d1] should be constrained to the dofs E and F at [e0, e1] and [f0, f1] as \(D = \alpha E + \beta F\) the dictionary should be:

{numpy.array([d0, d1], dtype=mesh.geometry.x.dtype).tobytes():
    {numpy.array([e0, e1], dtype=mesh.geometry.x.dtype).tobytes(): alpha,
    numpy.array([f0, f1], dtype=mesh.geometry.x.dtype).tobytes(): beta}}
create_periodic_constraint_geometrical(V: FunctionSpace, indicator: Callable[[ndarray], ndarray], relation: Callable[[ndarray], ndarray], bcs: List[DirichletBC], scale: float32 | float64 | complex128 | complex64 = 1.0)[source]#

Create a periodic condition for all degrees of freedom whose physical location satisfies \(indicator(x_i)==True\), i.e. \(u(x_i) = scale * u(relation(x_i))\) for all \(x_i\)

Parameters:
  • V – The function space to assign the condition to. Should either be the space of the MPC or a sub space.

  • indicator – Lambda-function to locate degrees of freedom that should be slaves

  • relation – Lambda-function describing the geometrical relation to master dofs

  • bcs – Dirichlet boundary conditions for the problem (Periodic constraints will be ignored for these dofs)

  • scale – Float for scaling bc

create_periodic_constraint_topological(V: FunctionSpace, meshtag: MeshTags, tag: int, relation: Callable[[ndarray], ndarray], bcs: List[DirichletBC], scale: float32 | float64 | complex128 | complex64 = 1.0)[source]#

Create periodic condition for all closure dofs of on all entities in meshtag with value tag. \(u(x_i) = scale * u(relation(x_i))\) for all of \(x_i\) on marked entities.

Parameters:
  • V – The function space to assign the condition to. Should either be the space of the MPC or a sub space. meshtag: MeshTag for entity to apply the periodic condition on

  • tag – Tag indicating which entities should be slaves

  • relation – Lambda-function describing the geometrical relation

  • bcs – Dirichlet boundary conditions for the problem (Periodic constraints will be ignored for these dofs)

  • scale – Float for scaling bc

create_slip_constraint(space: FunctionSpace, facet_marker: Tuple[MeshTags, int], v: Function, bcs: List[DirichletBC] = [])[source]#

Create a slip constraint \(u \cdot v=0\) over the entities defined in facet_marker with the given index.

Parameters:
  • space – Function space (possible sub space) for the current constraint

  • facet_marker – Tuple containomg the mesh tag and marker used to locate degrees of freedom

  • v – Function containing the directional vector to dot your slip condition (most commonly a normal vector)

  • bcs – List of Dirichlet BCs (slip conditions will be ignored on these dofs)

Examples

Create constaint \(u\cdot n=0\) of all indices in mt marked with i

V = dolfinx.fem.functionspace(mesh, ("CG", 1))
mpc = MultiPointConstraint(V)
n = dolfinx.fem.Function(V)
mpc.create_slip_constaint(V, (mt, i), n)

Create slip constaint for a mixed function space:

cellname = mesh.ufl_cell().cellname()
Ve = basix.ufl.element(basix.ElementFamily.P, cellname , 2, shape=(mesh.geometry.dim,))
Qe = basix.ufl.element(basix.ElementFamily.P, cellname , 1)
me = basix.ufl.mixed_element([Ve, Qe])
W = dolfinx.fem.functionspace(mesh, me)
mpc = MultiPointConstraint(W)
n_space, _ = W.sub(0).collapse()
normal = dolfinx.fem.Function(n_space)
mpc.create_slip_constraint(W.sub(0), (mt, i), normal, bcs=[])

A slip condition cannot be applied on the same degrees of freedom as a Dirichlet BC, and therefore any Dirichlet bc for the space of the multi point constraint should be supplied.

cellname = mesh.ufl_cell().cellname()
Ve = basix.ufl.element(basix.ElementFamily.P, cellname , 2, shape=(mesh.geometry.dim,))
Qe = basix.ufl.element(basix.ElementFamily.P, cellname , 1)
me = basix.ufl.mixed_element([Ve, Qe])
W = dolfinx.fem.functionspace(mesh, me)
mpc = MultiPointConstraint(W)
n_space, _ = W.sub(0).collapse()
normal = Function(n_space)
bc = dolfinx.fem.dirichletbc(inlet_velocity, dofs, W.sub(0))
mpc.create_slip_constraint(W.sub(0), (mt, i), normal, bcs=[bc])
finalize() None[source]#

Finializes the multi point constraint. After this function is called, no new constraints can be added to the constraint. This function creates a map from the cells (local to index) to the slave degrees of freedom and builds a new index map and function space where unghosted master dofs are added as ghosts.

property function_space#

Return the function space for the multi-point constraint with the updated index map

homogenize(u: Function) None[source]#

For a vector, homogenize (set to zero) the vector components at the multi-point constraint slave DoF indices. This is particularly useful for nonlinear problems.

Parameters:

u – The input vector

property is_slave: ndarray#

Returns a vector of integers where the ith entry indicates if a degree of freedom (local to process) is a slave.

property masters: AdjacencyList_int32#

Returns an adjacency-list whose ith node corresponds to a degree of freedom (local to process), and links the corresponding master dofs (local to process).

Examples

masters = mpc.masters
masters_of_dof_i = masters.links(i)
property num_local_slaves#

Return the number of slaves owned by the current process.

property slaves#

Returns the degrees of freedom for all slaves local to process

dolfinx_mpc.apply_lifting(b: Vec, form: List[Form], bcs: List[List[DirichletBC]], constraint: MultiPointConstraint, x0: List[Vec] = [], scale: float32 | float64 | complex128 | complex64 = 1.0)[source]#

Apply lifting to vector b, i.e. \(b = b - scale \cdot K^T (A_j (g_j - x0_j))\)

Parameters:
  • b – PETSc vector to assemble into

  • form – The linear form

  • bcs – List of Dirichlet boundary conditions

  • constraint – The multi point constraint

  • x0 – List of vectors

  • scale – Scaling for lifting

dolfinx_mpc.assemble_matrix(form: Form, constraint: MultiPointConstraint | Sequence[MultiPointConstraint], bcs: Sequence[DirichletBC] | None = None, diagval: float64 = 1, A: Mat | None = None) Mat[source]#

Assemble a compiled DOLFINx bilinear form into a PETSc matrix with corresponding multi point constraints and Dirichlet boundary conditions.

Parameters:
  • form – The compiled bilinear variational form

  • constraint – The multi point constraint

  • bcs – Sequence of Dirichlet boundary conditions

  • diagval – Value to set on the diagonal of the matrix

  • A – PETSc matrix to assemble into

Returns:

The matrix with the assembled bi-linear form #type: ignore

Return type:

_PETSc.Mat

dolfinx_mpc.assemble_matrix_nest(A: Mat, a: Sequence[Sequence[Form]], constraints: Sequence[MultiPointConstraint], bcs: Sequence[DirichletBC] = [], diagval: float64 = 1)[source]#

Assemble a compiled DOLFINx bilinear form into a PETSc matrix of type “nest” with corresponding multi point constraints and Dirichlet boundary conditions.

Parameters:
  • a – The compiled bilinear variational form provided in a rank 2 list

  • constraints – An ordered list of multi point constraints

  • bcs – Sequence of Dirichlet boundary conditions

  • diagval – Value to set on the diagonal of the matrix (Default 1)

  • A – PETSc matrix to assemble into

dolfinx_mpc.assemble_vector(form: Form, constraint: MultiPointConstraint, b: Vec | None = None) Vec[source]#

Assemble a linear form into vector b with corresponding multi point constraint

Parameters:
  • form – The linear form

  • constraint – The multi point constraint

  • b – PETSc vector to assemble

Returns:

The vector with the assembled linear form (b if supplied)

dolfinx_mpc.assemble_vector_nest(b: Vec, L: Sequence[Form], constraints: Sequence[MultiPointConstraint])[source]#

Assemble a linear form into a PETSc vector of type “nest”

Parameters:
  • b – A PETSc vector of type “nest”

  • L – A sequence of linear forms

  • constraints – An ordered list of multi point constraints

dolfinx_mpc.create_matrix_nest(a: Sequence[Sequence[Form]], constraints: Sequence[MultiPointConstraint])[source]#

Create a PETSc matrix of type “nest” with appropriate sparsity pattern given the provided multi points constraints

Parameters:

a – The compiled bilinear variational form provided in a rank 2 list constraints: An ordered list of multi point constraints

dolfinx_mpc.create_sparsity_pattern(form: Form, mpc: MultiPointConstraint | Sequence[MultiPointConstraint])[source]#

Create sparsity-pattern for MPC given a compiled DOLFINx form

Parameters:
  • form – The form

  • mpc – For square forms, the MPC. For rectangular forms a list of 2 MPCs on axis 0 & 1, respectively

dolfinx_mpc.create_vector_nest(L: Sequence[Form], constraints: Sequence[MultiPointConstraint]) Vec[source]#

Create a PETSc vector of type “nest” appropriate for the provided multi point constraints

Parameters:
  • L – A sequence of linear forms

  • constraints – An ordered list of multi point constraints

Returns:

A PETSc vector of type “nest” #type: ignore

Return type:

PETSc.Vec