-
-
Notifications
You must be signed in to change notification settings - Fork 30
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
Conflicting __eq__
/__hash__
in basix.ufl._BlockedElement
#717
Comments
A further note: the inconsistency between import basix, basix.ufl
import dolfinx as dfx
import ufl
from mpi4py import MPI
msh = dfx.mesh.create_unit_square(MPI.COMM_WORLD, 4, 4)
elem_scalar = basix.ufl.element(basix.ElementFamily.P, msh.basix_cell(), 1, gdim=2)
elem_vector = basix.ufl.blocked_element(elem_scalar, shape=(2,), gdim=2)
elem_tensor = basix.ufl.blocked_element(elem_scalar, shape=(2, 2), symmetry=True, gdim=2)
elem_mixed = basix.ufl.mixed_element([elem_vector, elem_tensor])
V = dfx.fem.functionspace(msh, elem_mixed)
u1, u2 = ufl.TrialFunctions(V)
v1, v2 = ufl.TestFunctions(V)
L = ufl.inner(ufl.grad(u1), v2) * ufl.dx
dfx.fem.form(L) producing
|
Yes, gdim should be added to the The other updates you suggest sound correct too. Do you want to open a PR to fix these? |
Sure, I'll have a go at it. |
I've opened a draft pull request with an initial implementation. Would it be good to add a test (maybe in |
Yes, it would be sensible to add a test there |
Describe the issue
Using the
gdim
argument inbasix.ufl.element
leads to conflicting results from_BlockedElement.__eq__
and_BlockedElement.__hash__
, as illustrated in the following MWE.results in
Expected behavior
This MWE should run successfully, since the hash values should be equal whenever
__eq__
returns True.Though not strictly required, I would expect them to be different when
__eq__
returns False.Proposed solution(s)
The solution could take one of the following two forms, depending on the intended functionality of the library:
_BlockedElement.__eq__
incorporating thegdim
. This is the correct solution if two elements identical except forgdim
should be considered different. (I.e.elem1 == elem2
evaluates toFalse
in the above MWE.) I expect this is the desired solution, but wanted to confirm that this is the intended functionality.gdim
in therepr
of_BlockedElement
, since the hash is based onrepr
. (I.e.elem1 == elem2
evaluates toTrue
in the above MWE.)We may also need to take a look at
_BasixElement
,_ComponentElement
, and_MixedElement
to make sure they are handlinggdim
correctly. Furthermore, we might need to check_QuadratureElement
for a similar issue, since the cell name and the pullback are included in therepr
but not in__eq__
.Also related, the
repr
of_RealElement
is wrong. Hopefully this MWE illustrates it sufficiently 😂 :producing:
The text was updated successfully, but these errors were encountered: