Skip to content

Refactor dolfinx::fem::Form to use local indexing of forms, rather than integral ids #3740

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Open
wants to merge 10 commits into
base: main
Choose a base branch
from

Conversation

jorgensd
Copy link
Member

@jorgensd jorgensd commented May 23, 2025

As shown in: #3735, there are many UFL forms (especially in the case of mixed quadrature rules), where the subdomain_id is not mapped 1-1 with an integration kernel.

This PR refactors the Form class, to rather use the local indexing from the generated code (where integrals have been grouped) as the lookup key for the integral.

Only a minor change occurs in the user-interface, as Form::subdomain_ids no longer maps to the UFL subdomain_ids.
However, this is currently only used for internal looping, as the integral_idx -> subdomain_id mapping happens within the create_form_factory function.

This allows for different quadrature degrees in different parts of a form.

Example:

import basix.ufl
import ufl

cell = basix.CellType.triangle
coord_element = basix.ufl.element("Lagrange", cell, 1, shape=(2,))
domain = ufl.Mesh(coord_element)

degree = 2
el = basix.ufl.element("Lagrange", cell, degree)
V = ufl.FunctionSpace(domain, el)
u = ufl.TrialFunction(V)
v = ufl.TestFunction(V)

dx = ufl.Measure("dx", domain=domain)

a_mass = ufl.inner(u, v)*dx((1, 2), degree=2*degree)
a_stiffness = ufl.inner(ufl.grad(u), ufl.grad(v))*dx((1,), degree=2*(degree-1))

a = a_mass + a_stiffness

This generates two different integration kernels in FFCx, with the integral information as:

uint64_t finite_element_hashes_form_e26395feca75c800a7470ee9a844c8cdc67aff23[2] = {UINT64_C(16933917890882727401), UINT64_C(16933917890882727401)};
int form_integral_offsets_form_e26395feca75c800a7470ee9a844c8cdc67aff23[4] = {0, 3, 3, 3};
static ufcx_integral* form_integrals_form_e26395feca75c800a7470ee9a844c8cdc67aff23[3] = {&integral_0dcf3cd8d98806a1645751b277f55c83190bebb9_triangle, &integral_be233757c48f146683317dc29c6fe8fa75e720a8_triangle, &integral_be233757c48f146683317dc29c6fe8fa75e720a8_triangle};
int form_integral_ids_form_e26395feca75c800a7470ee9a844c8cdc67aff23[3] = {1, 1, 2};

As seen in the latter line here, form_integral_ids_form_* is not a list of unique integers.
In the main branch of DOLFINx, we assume that this id is unique, which is the root cause of #3735.

@jorgensd jorgensd requested a review from chrisrichardson May 26, 2025 07:17
@garth-wells
Copy link
Member

Can you add some gentle background to the PR description?

@jorgensd
Copy link
Member Author

Can you add some gentle background to the PR description?

I've added an example to the description of the PR.

@jorgensd jorgensd added the enhancement New feature or request label May 29, 2025
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
enhancement New feature or request high-priority
Projects
None yet
Development

Successfully merging this pull request may close these issues.

[BUG]: Problem when assembling different integrals for different tags
2 participants