qugar.dolfinx

Custom DOLFINx forms for runtime quadratures

Functions

form_custom(form[, dtype, ...])

Creates a CustomForm or a list of CustomForm to be integrated over custom domains.

mapped_normal(domain[, normalize])

Returns a normal vector of a custom unfitted boundary mapped with the domain's geometry.

Classes

CustomForm(form, domain, itg_data[, ...])

Form for custom integrals.

ds_bdry_unf(domain[, subdomain_id, ...])

This is a new ufl Measure class for an unfitted custom boundary.

LinearProblem(a, L[, bcs, u, petsc_options, ...])

Class for solving a linear variational problem.

NonlinearProblem(F, u[, bcs, J, ...])

Nonlinear problem class for solving the non-linear problems.

class qugar.dolfinx.CustomForm(form: Form_complex64 | Form_complex128 | Form_float32 | Form_float64, domain: UnfittedDomainABC, itg_data: list[IntegralData], ufcx_form=None, code: str | None = None, module: ModuleType | None = None)[source]

Bases: Form

Form for custom integrals.

It derives from dolfinx.forms.Form just adding an extra functionality for computing the required custom coefficients at runtime (the method pack_coefficients).

A custom finite element form.

Note

CustomForms should normally be constructed using form_custom() and not using this class initialiser. This class is combined with different base classes that depend on the scalar type used in the Form.

Parameters:
  • form – Compiled form object.

  • ufcx_form – UFCx form.

  • code – Form C++ code.

  • module – CFFI module.

pack_coefficients() dict[tuple[IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]][source]

Function for generating the coefficients needed for computing custom integrals at runtime.

Note

This function mimics the behaviour of dolfinx.cpp.fem.pack_coefficients.

Returns:

Generated custom coefficients.

Return type:

dict[tuple[IntegralType, int], npt.NDArray]

update_coefficients(old_coeffs: dict[tuple[IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]]) dict[tuple[IntegralType, int], ndarray[Any, dtype[_ScalarType_co]]][source]

Function for updating the coefficients needed for computing custom integrals at runtime.

It updated coefficients previously generated by just recomputing the part of the coefficients associated to the DOLFINx coefficients and keeping the part of the coefficients associated to the custom integrals.

Returns:

Generated custom coefficients.

Return type:

dict[tuple[IntegralType, int], npt.NDArray]

class qugar.dolfinx.LinearProblem(a: Form, L: Form, bcs: list[DirichletBC] = [], u: Function | None = None, petsc_options: dict | None = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]

Bases: LinearProblem

Class for solving a linear variational problem.

Analogous to DOLFINx’ LinearProblem, but for unfitted domains.

Solves of the form \(a(u, v) = L(v) \, \forall v \in V\) using PETSc as a linear algebra backend.

Initialize solver for a linear variational problem defined over an unfitted domain.

Parameters:
  • a (ulf.Form) – A bilinear UFL form, the left hand side of the variational problem. This form must be defined over an uniftted domain.

  • L (ulf.Form) – A linear UFL form, the right hand side of the variational problem. This form must be defined over an uniftted domain.

  • bcs (list[DirichletBC]) – A list of Dirichlet boundary conditions.

  • u (Optional[_Function]) – The solution function. It will be created if not provided.

  • petsc_options (Optional[dict]) – Options that are passed to the linear algebra backend PETSc. For available choices for the ‘petsc_options’ kwarg, see the PETSc documentation.

  • form_compiler_options (Optional[dict]) – Options used in FFCx compilation of this form. Run ffcx --help at the commandline to see all available options.

  • jit_options (Optional[dict]) – Options used in CFFI JIT compilation of C code generated by FFCx See DOLFINx’s python/dolfinx/jit.py for all available options. Takes priority over all other option values.

Example:

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

Solve the problem.

class qugar.dolfinx.NonlinearProblem(F: Form, u: Function, bcs: list[DirichletBC] = [], J: Form = None, form_compiler_options: dict | None = None, jit_options: dict | None = None)[source]

Bases: NonlinearProblem

Nonlinear problem class for solving the non-linear problems.

Solves problems of the form \(F(u, v) = 0 \ \forall v \in V\) using PETSc as the linear algebra backend.

Initialize solver for solving a non-linear problem using Newton’s method`.

Parameters:
  • F – The PDE residual F(u, v)

  • u – The unknown

  • bcs – List of Dirichlet boundary conditions

  • J – UFL representation of the Jacobian (Optional)

  • form_compiler_options – Options used in FFCx compilation of this form. 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 python/dolfinx/jit.py for all available options. Takes priority over all other option values.

Example:

problem = LinearProblem(F, u, [bc0, bc1])
F(x: Vec, b: Vec) None[source]

Assemble the residual F into the vector b.

Parameters:
  • x – The vector containing the latest solution

  • b – Vector to assemble the residual into

J(x: Vec, A: Mat) None[source]

Assemble the Jacobian matrix.

Parameters:

x – The vector containing the latest solution

class qugar.dolfinx.ds_bdry_unf(domain: AbstractDomain | Mesh, subdomain_id: str | int | tuple[int] = 'everywhere', metadata: dict | None = None, subdomain_data: MeshTags | list[tuple[int, ndarray[Any, dtype[int32]]]] | None = None, degree: int | None = None)[source]

Bases: Measure

This is a new ufl Measure class for an unfitted custom boundary.

It has the same functionalities as ufl.dx with the only difference that when multiplied by an integrand, it introduces the necessary correction for accounting for the boundary orientation by using its normal.

In order to differentiate the generated measure from others, sets the option unfitted_custom_boundary equal to True in the quadrature’s metadata, and uses a custom quadrature with two (fake) points.

Initialize.

Parameters:
  • domain (ufl.AbstractDomain | dolfinx.mesh.Mesh) – An AbstractDomain object (most often a Mesh).

  • subdomain_id (str | int | tuple[int], optional) – either string “everywhere”, a single subdomain id int, or tuple of ints. Defaults to “everywhere”.

  • metadata (dict | None, optional) – Dictionary with additional compiler-specific parameters for optimization or debugging of generated code. Defaults to None.

  • subdomain_data (dolfinx.mesh.MeshTags | list[tuple[int, npt.NDArray[np.int32]]] | None, optional) – Object representing data to interpret subdomain_id with. Defaults to None.

  • degree (int, optional) – The degree of the quadrature rule.

qugar.dolfinx.form_custom(form: typing.Union[ufl.Form, typing.Iterable[ufl.Form]], dtype: npt.DTypeLike = <class 'numpy.float64'>, form_compiler_options: typing.Optional[dict] = None, jit_options: typing.Optional[dict] = None, entity_maps: typing.Optional[dict[Mesh, npt.NDArray[np.int32]]] = None) CustomForm | None | list[CustomForm | None][source]

Creates a CustomForm or a list of CustomForm to be integrated over custom domains.

Note

This function is just a copy of dolfinx.fem.forms.form with some small modifications. Namely:

  • replacing dolfinx.jit.ffcx_jit with qugar.dolfinx.jit.ffcx_jit.

  • creating a new form class CustomForm (that derives from dolfinx.forms.Form and adds a new functionalitiy for computing the required custom coefficients at runtime).

Parameters:
  • form – A UFL form or list(s) of UFL forms.

  • domain (UnfittedDomainABC) – The unfitted domain to integrate over.

  • dtype – Scalar type to use for the compiled form.

  • form_compiler_options – See ffcx_jit

  • jit_options – See ffcx_jit.

  • entity_maps – If any trial functions, test functions, or coefficients in the form are not defined over the same mesh as the integration domain, entity_maps must be supplied. For each key (a mesh, different to the integration domain mesh) a map should be provided relating the entities in the integration domain mesh to the entities in the key mesh e.g. for a key-value pair (msh, emap) in entity_maps, emap[i] is the entity in msh corresponding to entity i in the integration domain mesh.

Returns:

Compiled finite element CustomForm. It can be a single form or a list of forms if the input form is a list of forms. If the provided form is empty, it returns None.

Return type:

CustomForm | list[CustomForm]

Note

This function is responsible for the compilation of a UFL form (using FFCx) and attaching coefficients and domains specific data to the underlying C++ form. It dynamically create a Form instance with an appropriate base class for the scalar type, e.g. _cpp.fem.Form_float64().

qugar.dolfinx.mapped_normal(domain: AbstractDomain | Mesh, normalize: bool = True)[source]

Returns a normal vector of a custom unfitted boundary mapped with the domain’s geometry. I.e., the normal in the physical space.

Parameters:
  • domain (ufl.AbstractDomain | dolfinx.mesh.Mesh) – An AbstractDomain object (most often a Mesh) to which the normal is associated to.

  • normalize (bool) – If True, the returned normal is normalized, otherwise (if False) it is not normalized and its norm is the ratio of the deformed to reference area cells.

Returns:

Mapped normal.