qugar.dolfinx
Custom DOLFINx forms for runtime quadratures
Functions
|
Creates a CustomForm or a list of CustomForm to be integrated over custom domains. |
|
Returns a normal vector of a custom unfitted boundary mapped with the domain's geometry. |
Classes
|
Form for custom integrals. |
|
This is a new ufl Measure class for an unfitted custom boundary. |
- class qugar.dolfinx.CustomForm(form: Form_complex64 | Form_complex128 | Form_float32 | Form_float64, 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(quad_gen: QuadGenerator) dict[tuple[IntegralType, int], ndarray[Any, dtype[float32 | float64 | complex64 | complex128]]] [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
.- Parameters:
quad_gen (QuadGenerator) – Custom quadrature generator.
- Returns:
dict[ tuple[IntegralType, int], npt.NDArray[np.float32 | np.float64 | np.complex64 | np.complex128]]: Generated custom coefficients.
- class qugar.dolfinx.dx_bdry_unf(domain: AbstractDomain | Mesh, subdomain_id: str | int | tuple[int] = 'everywhere', metadata: dict | None = None, subdomain_data: MeshTags | 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 | 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)[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
withqugar.dolfinx.jit.ffcx_jit
.creating a new form class
CustomForm
(that derives fromdolfinx.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.
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.
- 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 (ifFalse
) it is not normalized and its norm is the ratio of the deformed to reference area cells.
- Returns:
Mapped normal.