dolfinx.fem.problems#

High-level problem classes using native linear algebra objects.

Users with advanced requirements should use dolfinx.fem.petsc.

Classes

LinearProblem(a, L, bcs, u, dtype, ...)

High-level class for solving linear problems using SuperLU_DIST.

class dolfinx.fem.problems.LinearProblem(a: Form, L: Form, bcs: Sequence[DirichletBC] | None = None, u: Function | None = None, dtype: DTypeLike = <class 'numpy.float64'>, superlu_dist_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None, entity_maps: Sequence[EntityMap] | None = None)[source]#

Bases: object

High-level class for solving linear problems using SuperLU_DIST.

Solves problems of the form \(a(u, v) = L(v) \; \forall v \in V\) using dolfinx.la.superlu_dist.SuperLUDistSolver as the linear solver.

Note

DOLFINx must be built with SuperLU_DIST to use this class.

Initialize SuperLU_DIST solver for a linear variational problem.

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

  • L – Linear UFL form, the right-hand side of the variational problem.

  • bcs – Dirichlet boundary conditions to apply to the variational problem.

  • u – Solution function. Created if not provided.

  • dtype – Scalar type for form compilation. Must match u.dtype if u is provided.

  • superlu_dist_options – Options passed to the SuperLU_DIST solver. See the SuperLU_DIST User’s Guide for available options and values.

  • form_compiler_options – Options used in FFCx compilation of all forms. Run ffcx --help at the command line to see all available options.

  • jit_options – Options used in CFFI JIT compilation of C code generated by FFCx. See dolfinx.jit.ffcx_jit() for all available options. Takes priority over all other option values.

  • entity_maps – If any trial functions, test functions, or coefficients in the form are not defined over the same mesh as the integration domain, a corresponding EntityMap must be provided.

Example:

problem = LinearProblem(a, L, bcs=bc,
    superlu_dist_options={"SymmetricMode": "YES"})
u_h = problem.solve()
property A: MatrixCSR#

Left-hand side matrix.

property L: Form#

The compiled linear form representing the right-hand side.

property a: Form#

The compiled bilinear form representing the left-hand side.

property b: Vector#

Right-hand side vector.

solve() Function[source]#

Solve the problem.

This method updates the solution function u stored in the problem instance.

Returns:

The solution function.

property u: Function#

Solution function.

Note

The function does not share memory with the solution vector x.

property x: Vector#

Solution vector.

Note

The vector does not share memory with the solution function u.