Skip to content

Conflicting __eq__/__hash__ in basix.ufl._BlockedElement #717

@conpierce8

Description

@conpierce8

Describe the issue

Using the gdim argument in basix.ufl.element leads to conflicting results from _BlockedElement.__eq__ and _BlockedElement.__hash__, as illustrated in the following MWE.

import basix, basix.ufl
elem1 = basix.ufl.element(basix.ElementFamily.P, basix.CellType.triangle, 1, shape=(2,))
elem2 = basix.ufl.element(basix.ElementFamily.P, basix.CellType.triangle, 1, shape=(2,), gdim=2)
assert (elem1 == elem2 and hash(elem1) == hash(elem2)) or (elem1 != elem2 and hash(elem1) != hash(elem2))

results in

Traceback (most recent call last):
  File "/root/mwe.py", line 4, in <module>
    assert (elem1 == elem2 and hash(elem1) == hash(elem2)) or (elem1 != elem2 and hash(elem1) != hash(elem2))
AssertionError

Expected behavior

This MWE should run successfully, since the hash values should be equal whenever __eq__ returns True.

The only required property is that objects which compare equal have the same hash value

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:

  1. Add an extra condition to _BlockedElement.__eq__ incorporating the gdim. This is the correct solution if two elements identical except for gdim should be considered different. (I.e. elem1 == elem2 evaluates to False in the above MWE.) I expect this is the desired solution, but wanted to confirm that this is the intended functionality.
        def __eq__(self, other) -> bool:
         """Check if two elements are equal."""
         return (
             isinstance(other, _BlockedElement) and self._block_size == other._block_size
             and self.block_shape == other.block_shape and self.sub_element == other.sub_element
             and self._gdim == other._gdim)
  2. Remove any reference to gdim in the repr of _BlockedElement, since the hash is based on repr. (I.e. elem1 == elem2 evaluates to True in the above MWE.)

We may also need to take a look at _BasixElement, _ComponentElement, and _MixedElement to make sure they are handling gdim correctly. Furthermore, we might need to check _QuadratureElement for a similar issue, since the cell name and the pullback are included in the repr but not in __eq__.

Also related, the repr of _RealElement is wrong. Hopefully this MWE illustrates it sufficiently 😂 :

import basix, basix.ufl
real_element = basix.ufl.real_element(basix.CellType.triangle, (2,))
repr(real_element)

producing:

'RealElement(<functools._lru_cache_wrapper object at 0x7f51cfb3e770>)'

Metadata

Metadata

Assignees

No one assigned

    Labels

    No labels
    No labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions