Skip to content
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

Add sobolev_space method to base finite element class #118

Merged
merged 10 commits into from
Aug 16, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension


Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 4 additions & 4 deletions .github/workflows/fenicsx-tests.yml
Original file line number Diff line number Diff line change
Expand Up @@ -33,14 +33,14 @@ jobs:

- name: Install Basix
run: |
python3 -m pip install git+https://github.com/FEniCS/basix.git
python3 -m pip install git+https://github.com/FEniCS/basix.git@mscroggs/ufl-update

- name: Clone FFCx
uses: actions/checkout@v2
with:
path: ./ffcx
repository: FEniCS/ffcx
ref: main
ref: mscroggs/ufl_update
- name: Install FFCx
run: |
cd ffcx
Expand Down Expand Up @@ -78,8 +78,8 @@ jobs:

- name: Install Basix and FFCx
run: |
python3 -m pip install git+https://github.com/FEniCS/basix.git
python3 -m pip install git+https://github.com/FEniCS/ffcx.git
python3 -m pip install git+https://github.com/FEniCS/basix.git@mscroggs/ufl-update
python3 -m pip install git+https://github.com/FEniCS/ffcx.git@mscroggs/ufl_update

- name: Clone DOLFINx
uses: actions/checkout@v2
Expand Down
2 changes: 2 additions & 0 deletions .gitignore
Original file line number Diff line number Diff line change
Expand Up @@ -8,3 +8,5 @@ test/.pytest_cache
*.egg-info
dist/
release/

.*.swp
10 changes: 6 additions & 4 deletions ufl/finiteelement/brokenelement.py
Original file line number Diff line number Diff line change
Expand Up @@ -14,7 +14,6 @@ class BrokenElement(FiniteElementBase):
"""The discontinuous version of an existing Finite Element space."""
def __init__(self, element):
self._element = element
self._repr = "BrokenElement(%s)" % repr(element)

family = "BrokenElement"
cell = element.cell()
Expand All @@ -25,15 +24,18 @@ def __init__(self, element):
FiniteElementBase.__init__(self, family, cell, degree,
quad_scheme, value_shape, reference_value_shape)

def __repr__(self):
return f"BrokenElement({repr(self._element)})"

def mapping(self):
return self._element.mapping()

def reconstruct(self, **kwargs):
return BrokenElement(self._element.reconstruct(**kwargs))

def __str__(self):
return "BrokenElement(%s)" % str(self._element)
return f"BrokenElement({repr(self._element)})"

def shortstr(self):
"Format as string for pretty printing."
return "BrokenElement(%s)" % str(self._element.shortstr())
"""Format as string for pretty printing."""
return f"BrokenElement({repr(self._element)})"
11 changes: 7 additions & 4 deletions ufl/finiteelement/enrichedelement.py
Original file line number Diff line number Diff line change
Expand Up @@ -61,14 +61,11 @@ def __init__(self, *elements):
quad_scheme, value_shape,
reference_value_shape)

# Cache repr string
self._repr = "%s(%s)" % (class_name, ", ".join(repr(e) for e in self._elements))

def mapping(self):
return self._elements[0].mapping()

def sobolev_space(self):
"Return the underlying Sobolev space."
"""Return the underlying Sobolev space."""
elements = [e for e in self._elements]
if all(e.sobolev_space() == elements[0].sobolev_space()
for e in elements):
Expand Down Expand Up @@ -104,6 +101,9 @@ def is_cellwise_constant(self):
element is spatially constant over each cell."""
return all(e.is_cellwise_constant() for e in self._elements)

def __repr__(self):
return "EnrichedElement(" + ", ".join(repr(e) for e in self._elements) + ")"

def __str__(self):
"Format as string for pretty printing."
return "<%s>" % " + ".join(str(e) for e in self._elements)
Expand All @@ -127,6 +127,9 @@ def is_cellwise_constant(self):
element is spatially constant over each cell."""
return False

def __repr__(self):
return "NodalEnrichedElement(" + ", ".join(repr(e) for e in self._elements) + ")"

def __str__(self):
"Format as string for pretty printing."
return "<Nodal enriched element(%s)>" % ", ".join(str(e) for e in self._elements)
Expand Down
13 changes: 8 additions & 5 deletions ufl/finiteelement/finiteelement.py
Original file line number Diff line number Diff line change
Expand Up @@ -24,10 +24,8 @@
class FiniteElement(FiniteElementBase):
"The basic finite element class for all simple finite elements."
# TODO: Move these to base?
__slots__ = ("_short_name",
"_sobolev_space",
"_mapping",
"_variant")
__slots__ = ("_short_name", "_sobolev_space",
"_mapping", "_variant", "_repr")

def __new__(cls,
family,
Expand Down Expand Up @@ -188,11 +186,16 @@ def __init__(self,
repr(self.family()), repr(self.cell()), repr(self.degree()), quad_str, var_str)
assert '"' not in self._repr

def __repr__(self):
"""Format as string for evaluation as Python object."""
return self._repr

def mapping(self):
"""Return the mapping type for this element ."""
return self._mapping

def sobolev_space(self):
"Return the underlying Sobolev space."
"""Return the underlying Sobolev space."""
return self._sobolev_space

def variant(self):
Expand Down
35 changes: 17 additions & 18 deletions ufl/finiteelement/finiteelementbase.py
Original file line number Diff line number Diff line change
Expand Up @@ -15,24 +15,19 @@
from ufl.utils.dicts import EmptyDict
from ufl.log import error
from ufl.cell import AbstractCell, as_cell
from abc import ABC, abstractmethod


class FiniteElementBase(object):
class FiniteElementBase(ABC):
"Base class for all finite elements."
__slots__ = ("_family",
"_cell",
"_degree",
"_quad_scheme",
"_value_shape",
"_reference_value_shape",
"_repr",
"__weakref__")
__slots__ = ("_family", "_cell", "_degree", "_quad_scheme",
"_value_shape", "_reference_value_shape", "__weakref__")

# TODO: Not all these should be in the base class! In particular
# family, degree, and quad_scheme do not belong here.
def __init__(self, family, cell, degree, quad_scheme, value_shape,
reference_value_shape):
"Initialize basic finite element data."
"""Initialize basic finite element data."""
if not isinstance(family, str):
error("Invalid family type.")
if not (degree is None or isinstance(degree, (int, tuple))):
Expand All @@ -54,12 +49,20 @@ def __init__(self, family, cell, degree, quad_scheme, value_shape,
self._reference_value_shape = reference_value_shape
self._quad_scheme = quad_scheme

@abstractmethod
def __repr__(self):
"""Format as string for evaluation as Python object.
"""Format as string for evaluation as Python object."""
pass

NB! Assumes subclass has assigned its repr string to self._repr.
"""
return self._repr
@abstractmethod
def sobolev_space(self):
"""Return the underlying Sobolev space."""
pass

@abstractmethod
def mapping(self):
"""Return the mapping type for this element."""
pass

def _ufl_hash_data_(self):
return repr(self)
Expand Down Expand Up @@ -101,10 +104,6 @@ def quadrature_scheme(self):
"Return quadrature scheme of finite element."
return self._quad_scheme

def mapping(self):
"Not implemented."
error("Missing implementation of mapping().")

def cell(self):
"Return cell of finite element."
return self._cell
Expand Down
30 changes: 18 additions & 12 deletions ufl/finiteelement/hdivcurl.py
Original file line number Diff line number Diff line change
Expand Up @@ -14,11 +14,10 @@
class HDivElement(FiniteElementBase):
"""A div-conforming version of an outer product element, assuming
this makes mathematical sense."""
__slots__ = ("_element",)
__slots__ = ("_element", )

def __init__(self, element):
self._element = element
self._repr = "HDivElement(%s)" % repr(element)

family = "TensorProductElement"
cell = element.cell()
Expand All @@ -31,22 +30,25 @@ def __init__(self, element):
FiniteElementBase.__init__(self, family, cell, degree,
quad_scheme, value_shape, reference_value_shape)

def __repr__(self):
return f"HDivElement({repr(self._element)})"

def mapping(self):
return "contravariant Piola"

def sobolev_space(self):
"Return the underlying Sobolev space."
"""Return the underlying Sobolev space."""
return HDiv

def reconstruct(self, **kwargs):
return HDivElement(self._element.reconstruct(**kwargs))

def __str__(self):
return "HDivElement(%s)" % str(self._element)
return f"HDivElement({repr(self._element)})"

def shortstr(self):
"Format as string for pretty printing."
return "HDivElement(%s)" % str(self._element.shortstr())
"""Format as string for pretty printing."""
return f"HDivElement({self._element.shortstr()})"


class HCurlElement(FiniteElementBase):
Expand All @@ -56,7 +58,6 @@ class HCurlElement(FiniteElementBase):

def __init__(self, element):
self._element = element
self._repr = "HCurlElement(%s)" % repr(element)

family = "TensorProductElement"
cell = element.cell()
Expand All @@ -70,6 +71,9 @@ def __init__(self, element):
FiniteElementBase.__init__(self, family, cell, degree, quad_scheme,
value_shape, reference_value_shape)

def __repr__(self):
return f"HCurlElement({repr(self._element)})"

def mapping(self):
return "covariant Piola"

Expand All @@ -81,11 +85,11 @@ def reconstruct(self, **kwargs):
return HCurlElement(self._element.reconstruct(**kwargs))

def __str__(self):
return "HCurlElement(%s)" % str(self._element)
return f"HCurlElement({repr(self._element)})"

def shortstr(self):
"Format as string for pretty printing."
return "HCurlElement(%s)" % str(self._element.shortstr())
return f"HCurlElement({self._element.shortstr()})"


class WithMapping(FiniteElementBase):
Expand All @@ -95,7 +99,6 @@ class WithMapping(FiniteElementBase):
remapped = WithMapping(E, "identity")
"""
def __init__(self, wrapee, mapping):
self._repr = "WithMapping(%s, %s)" % (repr(wrapee), mapping)
if mapping == "symmetries":
raise ValueError("Can't change mapping to 'symmetries'")
self._mapping = mapping
Expand All @@ -108,6 +111,9 @@ def __getattr__(self, attr):
raise AttributeError("'%s' object has no attribute '%s'" %
(type(self).__name__, attr))

def __repr__(self):
return f"WithMapping({repr(self.wrapee)}, {self._mapping})"

def value_shape(self):
gdim = self.cell().geometric_dimension()
mapping = self.mapping()
Expand Down Expand Up @@ -137,7 +143,7 @@ def reconstruct(self, **kwargs):
return type(self)(wrapee, mapping)

def __str__(self):
return "WithMapping(%s, mapping=%s)" % (self.wrapee, self._mapping)
return f"WithMapping({repr(self.wrapee)}, {self._mapping})"

def shortstr(self):
return "WithMapping(%s, %s)" % (self.wrapee.shortstr(), self._mapping)
return f"WithMapping({self.wrapee.shortstr()}, {self._mapping})"
26 changes: 17 additions & 9 deletions ufl/finiteelement/mixedelement.py
Original file line number Diff line number Diff line change
Expand Up @@ -90,10 +90,8 @@ def __init__(self, *elements, **kwargs):
FiniteElementBase.__init__(self, "Mixed", cell, degree, quad_scheme,
value_shape, reference_value_shape)

# Cache repr string
if type(self) is MixedElement:
self._repr = "MixedElement(%s)" % (
", ".join(repr(e) for e in self._sub_elements),)
def __repr__(self):
return "MixedElement(" + ", ".join(repr(e) for e in self._sub_elements) + ")"

def reconstruct_from_elements(self, *elements):
"Reconstruct a mixed element from new subelements."
Expand Down Expand Up @@ -125,6 +123,9 @@ def symmetry(self):
error("Size mismatch in symmetry algorithm.")
return sm or EmptyDict

def sobolev_space(self):
return max(e.sobolev_space() for e in self._sub_elements)

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Is this the right thing? Suppose I have CG1 x DG0. I think this will say that the whole thing is in H1

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Yes, looks like max might be what we actually want here?

Copy link
Collaborator

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think that's a safe conservative option. It might disallow implicit restrictions in some cases where they would actually work, but that can be dealt with in a later pass if someone cares enough.

def mapping(self):
if all(e.mapping() == "identity" for e in self._sub_elements):
return "identity"
Expand Down Expand Up @@ -247,6 +248,8 @@ def shortstr(self):
class VectorElement(MixedElement):
"A special case of a mixed finite element where all elements are equal."

__slots__ = ("_repr", "_mapping", "_sub_element")

def __init__(self, family, cell=None, degree=None, dim=None,
form_degree=None, quad_scheme=None, variant=None):
"""
Expand Down Expand Up @@ -310,8 +313,10 @@ def __init__(self, family, cell=None, degree=None, dim=None,
var_str = ", variant='" + variant + "'"

# Cache repr string
self._repr = "VectorElement(%s, dim=%d%s)" % (
repr(sub_element), len(self._sub_elements), var_str)
self._repr = f"VectorElement({repr(sub_element)}, dim={dim}{var_str})"

def __repr__(self):
return self._repr

def reconstruct(self, **kwargs):
sub_element = self._sub_element.reconstruct(**kwargs)
Expand Down Expand Up @@ -343,7 +348,7 @@ class TensorElement(MixedElement):
__slots__ = ("_sub_element", "_shape", "_symmetry",
"_sub_element_mapping",
"_flattened_sub_element_mapping",
"_mapping")
"_mapping", "_repr")

def __init__(self, family, cell=None, degree=None, shape=None,
symmetry=None, quad_scheme=None, variant=None):
Expand Down Expand Up @@ -447,8 +452,11 @@ def __init__(self, family, cell=None, degree=None, shape=None,
var_str = ", variant='" + variant + "'"

# Cache repr string
self._repr = "TensorElement(%s, shape=%s, symmetry=%s%s)" % (
repr(sub_element), repr(self._shape), repr(self._symmetry), var_str)
self._repr = (f"TensorElement({repr(sub_element)}, shape={shape}, "
f"symmetry={symmetry}{var_str})")

def __repr__(self):
return self._repr

def variant(self):
"""Return the variant used to initialise the element."""
Expand Down
8 changes: 6 additions & 2 deletions ufl/finiteelement/restrictedelement.py
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@

from ufl.finiteelement.finiteelementbase import FiniteElementBase
from ufl.log import error
from ufl.sobolevspace import L2

valid_restriction_domains = ("interior", "facet", "face", "edge", "vertex")

Expand All @@ -35,8 +36,11 @@ def __init__(self, element, restriction_domain):

self._restriction_domain = restriction_domain

self._repr = "RestrictedElement(%s, %s)" % (
repr(self._element), repr(self._restriction_domain))
def __repr__(self):
return f"RestrictedElement({repr(self._element)}, {repr(self._restriction_domain)})"

def sobolev_space(self):
return L2

def is_cellwise_constant(self):
"""Return whether the basis functions of this element is spatially
Expand Down
Loading