Skip to content

Conversation

@joewallwork
Copy link
Member

@joewallwork joewallwork commented Dec 5, 2025

Closes #3.

This PR moves over the norm and errornorm functions from Animate, which overload the Firedrake equivalents to account for vector l1, l2, and linf norms, as well as integral norms.

@joewallwork joewallwork marked this pull request as ready for review December 5, 2025 10:10
@joewallwork joewallwork moved this from In Progress to Ready for review in Mesh Adaptation development board Dec 5, 2025
if isinstance(v, fd.Cofunction):
v = cofunction2function(v)
condition = condition or fd.Constant(1.0)
norm_codes = {"l1": 0, "l2": 2, "linf": 3}
Copy link
Collaborator

Choose a reason for hiding this comment

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

Suggested change
norm_codes = {"l1": 0, "l2": 2, "linf": 3}
norm_codes = {
"l1": PETSc.NormType.NORM_1,
"l2": PETSc.NormType.NORM_2,
"linf": PETSc.NormType.NORM_INFINITY
}

v.interpolate(condition * v)
with v.dat.vec_ro as vv:
if norm_type == "Linf":
return vv.max()[1]
Copy link
Collaborator

Choose a reason for hiding this comment

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

Does this need an abs()? I presume "Linf" is in fact the same as "linf" ?

Comment on lines +73 to +78
integrand = {
"h1": lambda w: ufl.inner(w, w) + ufl.inner(ufl.grad(w), ufl.grad(w)),
"hdiv": lambda w: ufl.inner(w, w) + ufl.div(w) * ufl.div(w),
"hcurl": lambda w: ufl.inner(w, w)
+ ufl.inner(ufl.curl(w), ufl.curl(w)),
}[norm_type.lower()](v)
Copy link
Collaborator

Choose a reason for hiding this comment

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

I'm probably missing something, but what's the point of the lambda's? If you're worried about the cost of ufl symbolic assembly - I would just change it to an if block (or match case)

if boundary:
not_impl_err = "lp errors on the boundary not yet implemented."
raise NotImplementedError(not_impl_err)
v.interpolate(condition * v)
Copy link
Collaborator

Choose a reason for hiding this comment

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

Hm, this clobbers the input Function v - I'm not sure I like that...

# Case 1: point-wise norms
if norm_type[0] == "l":
v = u
v -= uh
Copy link
Collaborator

Choose a reason for hiding this comment

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

If I'm not mistaken -= for Functions is overloaded and basically does v.assign(v -uh) which again means that the values of the input u are overwritten as v is the same as u at this point. So then if you call this errornorm twice with the same u and uh, you get the wrong thing on the 2nd call.

:returns: the error norm value
:rtype: :class:`float`

Any other keyword arguments are passed to :func:`fd.norms.errornorm`.
Copy link
Collaborator

Choose a reason for hiding this comment

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

That's not true it seems

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

enhancement New feature or request

Projects

Status: Ready for review

Development

Successfully merging this pull request may close these issues.

Add norms module

3 participants