errornorm

dolfin.fem.norms.errornorm(u, uh, norm_type='l2', degree_rise=3, mesh=None)

Compute and return the error \(e = u - u_h\) in the given norm.

Arguments
u, uh
Functions
norm_type
Type of norm. The \(L^2\) -norm is default. For other norms, see norm.
degree_rise
The number of degrees above that of u_h used in the interpolation; i.e. the degree of piecewise polynomials used to approximate \(u\) and \(u_h\) will be the degree of \(u_h\) + degree_raise.
mesh
Optional Mesh on which to compute the error norm.

In simple cases, one may just define

e = u - uh

and evalute for example the square of the error in the \(L^2\) -norm by

assemble(e**2*dx(mesh))

However, this is not stable w.r.t. round-off errors considering that the form compiler may expand(#) the expression above to:

e**2*dx = u**2*dx - 2*u*uh*dx + uh**2*dx

and this might get further expanded into thousands of terms for higher order elements. Thus, the error will be evaluated by adding a large number of terms which should sum up to something close to zero (if the error is small).

This module computes the error by first interpolating both \(u\) and \(u_h\) to a common space (of high accuracy), then subtracting the two fields (which is easy since they are expressed in the same basis) and then evaluating the integral.

(#) If using the tensor representation optimizations. The quadrature represenation does not suffer from this problem.