Skip to content

Commit

Permalink
Merge pull request #28 from firedrakeproject/fabianl1908/variants
Browse files Browse the repository at this point in the history
Change default variant of Nedelec, BDM and RT elements
  • Loading branch information
wence- authored Mar 16, 2022
2 parents 5037645 + d1fc141 commit 5384beb
Show file tree
Hide file tree
Showing 6 changed files with 22 additions and 26 deletions.
10 changes: 5 additions & 5 deletions FIAT/brezzi_douglas_marini.py
Original file line number Diff line number Diff line change
Expand Up @@ -109,15 +109,15 @@ class BrezziDouglasMarini(finite_element.CiarletElement):
variant can be chosen from ["point", "integral", "integral(quadrature_degree)"]
"point" -> dofs are evaluated by point evaluation. Note that this variant has suboptimal
convergence order in the H(div)-norm
"integral" -> dofs are evaluated by quadrature rule. The quadrature degree is chosen to integrate
polynomials of degree 5*k so that most expressions will be interpolated exactly. This is important
when you want to have (nearly) divergence-preserving interpolation.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree
"integral" -> dofs are evaluated by quadrature rule of degree k.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree. You might
want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important
when you want to have (nearly) div-preserving interpolation.
"""

def __init__(self, ref_el, k, variant=None):

(variant, quad_deg) = check_format_variant(variant, k, "BDM")
(variant, quad_deg) = check_format_variant(variant, k)

if k < 1:
raise Exception("BDM_k elements only valid for k >= 1")
Expand Down
10 changes: 3 additions & 7 deletions FIAT/check_format_variant.py
Original file line number Diff line number Diff line change
@@ -1,19 +1,15 @@
import re
import warnings


def check_format_variant(variant, degree, element):
def check_format_variant(variant, degree):
if variant is None:
variant = "point"
warnings.simplefilter('always', DeprecationWarning)
warnings.warn('Variant of ' + element + ' element will change from point evaluation to integral evaluation.'
' You should project into variant="integral"', DeprecationWarning)
variant = "integral"

match = re.match(r"^integral(?:\((\d+)\))?$", variant)
if match:
variant = "integral"
quad_degree, = match.groups()
quad_degree = int(quad_degree) if quad_degree is not None else 5*(degree + 1)
quad_degree = int(quad_degree) if quad_degree is not None else (degree + 1)
if quad_degree < degree + 1:
raise ValueError("Warning, quadrature degree should be at least %s" % (degree + 1))
elif variant == "point":
Expand Down
8 changes: 4 additions & 4 deletions FIAT/nedelec.py
Original file line number Diff line number Diff line change
Expand Up @@ -360,17 +360,17 @@ class Nedelec(finite_element.CiarletElement):
variant can be chosen from ["point", "integral", "integral(quadrature_degree)"]
"point" -> dofs are evaluated by point evaluation. Note that this variant has suboptimal
convergence order in the H(curl)-norm
"integral" -> dofs are evaluated by quadrature rule. The quadrature degree is chosen to integrate
polynomials of degree 5*k so that most expressions will be interpolated exactly. This is important
"integral" -> dofs are evaluated by quadrature rule of degree k.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree. You might
want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important
when you want to have (nearly) curl-preserving interpolation.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree
"""

def __init__(self, ref_el, k, variant=None):

degree = k - 1

(variant, quad_deg) = check_format_variant(variant, degree, "Nedelec")
(variant, quad_deg) = check_format_variant(variant, degree)

if ref_el.get_spatial_dimension() == 3:
poly_set = NedelecSpace3D(ref_el, degree)
Expand Down
8 changes: 4 additions & 4 deletions FIAT/nedelec_second_kind.py
Original file line number Diff line number Diff line change
Expand Up @@ -227,15 +227,15 @@ class NedelecSecondKind(CiarletElement):
variant can be chosen from ["point", "integral", "integral(quadrature_degree)"]
"point" -> dofs are evaluated by point evaluation. Note that this variant has suboptimal
convergence order in the H(curl)-norm
"integral" -> dofs are evaluated by quadrature rule. The quadrature degree is chosen to integrate
polynomials of degree 5*k so that most expressions will be interpolated exactly. This is important
"integral" -> dofs are evaluated by quadrature rule of degree k.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree. You might
want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important
when you want to have (nearly) curl-preserving interpolation.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree
"""

def __init__(self, cell, k, variant=None):

(variant, quad_deg) = check_format_variant(variant, k, "Nedelec Second Kind")
(variant, quad_deg) = check_format_variant(variant, k)

# Check degree
assert k >= 1, "Second kind Nedelecs start at 1!"
Expand Down
10 changes: 5 additions & 5 deletions FIAT/raviart_thomas.py
Original file line number Diff line number Diff line change
Expand Up @@ -154,17 +154,17 @@ class RaviartThomas(finite_element.CiarletElement):
variant can be chosen from ["point", "integral", "integral(quadrature_degree)"]
"point" -> dofs are evaluated by point evaluation. Note that this variant has suboptimal
convergence order in the H(div)-norm
"integral" -> dofs are evaluated by quadrature rule. The quadrature degree is chosen to integrate
polynomials of degree 5*k so that most expressions will be interpolated exactly. This is important
when you want to have (nearly) divergence-preserving interpolation.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree
"integral" -> dofs are evaluated by quadrature rule of degree k.
"integral(quadrature_degree)" -> dofs are evaluated by quadrature rule of degree quadrature_degree. You might
want to choose a high quadrature degree to make sure that expressions will be interpolated exactly. This is important
when you want to have (nearly) div-preserving interpolation.
"""

def __init__(self, ref_el, k, variant=None):

degree = k - 1

(variant, quad_deg) = check_format_variant(variant, degree, "Raviart Thomas")
(variant, quad_deg) = check_format_variant(variant, degree)

poly_set = RTSpace(ref_el, degree)
dual = RTDualSet(ref_el, degree, variant, quad_deg)
Expand Down
2 changes: 1 addition & 1 deletion test/regression/fiat-reference-data-id
Original file line number Diff line number Diff line change
@@ -1 +1 @@
83d6c1d8f30d2c116398f496a4592ef541ea2843
0c8c97f7e4919402129e5ff3b54e3f0b9e902b7c

0 comments on commit 5384beb

Please sign in to comment.