1. Norms
Let \(f\) a bounded function on domain \(\Omega\).
1.1. L2 norms
Let \(f \in L^2(\Omega)\) you can evaluate the \(L^2\) norm using the normL2()
function:
1.1.1. Interface
normL2( _range, _expr, _quad, _geomap );
or squared norm:
normL2Squared( _range, _expr, _quad, _geomap );
Required parameters:
-
_range
= domain of integration -
_expr
= mesurable function
Optional parameters:
-
_quad
= quadrature to use.-
Default =
_Q<integer>()
-
-
_geomap
= type of geometric mapping.-
Default =
GEOMAP_OPT
-
1.1.2. Example
From doc/manual/laplacian/laplacian.cpp
double L2error =normL2( _range=elements( mesh ),
_expr=( idv( u )-g ) );
From doc/manual/stokes/stokes.cpp
mean
int main(int argc, char**argv )
{
Environment env( _argc=argc, _argv=argv,
_about=about(_name="mystokes",
_author="Feel++ Consortium",
_email="feelpp-devel@feelpp.org"));
// create the mesh
auto mesh = loadMesh(_mesh=new Mesh<Simplex< 2 > > );
// function space
auto Vh = THch<2>( mesh );
// element U=(u,p) in Vh
auto U = Vh->element();
auto u = U.element<0>();
auto p = U.element<1>();
// left hand side
auto a = form2( _trial=Vh, _test=Vh );
a = integrate(_range=elements(mesh),
_expr=trace(gradt(u)*trans(grad(u))) );
a+= integrate(_range=elements(mesh),
_expr=-div(u)*idt(p)-divt(u)*id(p));
auto syms = symbols<2>();
auto u1 = parse( option(_name="functions.alpha").as<std::string>(), syms );
auto u2 = parse( option(_name="functions.beta").as<std::string>(), syms );
matrix u_exact = matrix(2,1);
u_exact = u1,u2;
auto p_exact = parse( option(_name="functions.gamma").as<std::string>(), syms );
auto f = -laplacian( u_exact, syms ) + grad( p_exact, syms ).transpose();
LOG(INFO) << "rhs : " << f;
// right hand side
auto l = form1( _test=Vh );
l = integrate(_range=elements(mesh),
_expr=trans(expr<2,1,5>( f, syms ))*id(u));
a+=on(_range=boundaryfaces(mesh), _rhs=l, _element=u,
_expr=expr<2,1,5>(u_exact,syms));
// solve a(u,v)=l(v)
a.solve(_rhs=l,_solution=U);
double mean_p = mean(_range=elements(mesh),_expr=idv(p))(0,0);
double mean_p_exact = mean(_range=elements(mesh),_expr=expr(p_exact,syms))(0,0);
double l2error_u = normL2( _range=elements(mesh), _expr=idv(u)-expr<2,1,5>( u_exact, syms ) );
double l2error_p = normL2( _range=elements(mesh), _expr=idv(p)-mean_p-(expr( p_exact, syms )-mean_p_exact) );
LOG(INFO) << "L2 error norm u: " << l2error_u;
LOG(INFO) << "L2 error norm p: " << l2error_p;
// save results
auto e = exporter( _mesh=mesh );
e->add( "u", u );
e->add( "p", p );
e->save();
}
1.2. H^1 norm
In the same idea, you can evaluate the H1 norm or semi norm, for any function \(f \in H^1(\Omega)\):
where \(*\) is the scalar product \(\cdot\) when \(f\) is a scalar field and the frobenius scalar product \(:\) when \(f\) is a vector field.
1.2.1. Interface
normH1( _range, _expr, _grad_expr, _quad, _geomap );
or semi norm:
normSemiH1( _range, _grad_expr, _quad, _geomap );
Required parameters:
-
_range
= domain of integration -
_expr
= mesurable function -
_grad_expr
= gradient of function (Row vector!)
Optional parameters:
-
_quad
= quadrature to use.-
Default =
_Q<integer>()
-
-
_geomap
= type of geometric mapping.-
Default =
GEOMAP_OPT
-
normH1() returns a float containing the H^1 norm.
1.2.2. Example
With expression:
auto g = sin(2*pi*Px())*cos(2*pi*Py());
auto gradg = 2*pi*cos(2* pi*Px())*cos(2*pi*Py())*oneX()
-2*pi*sin(2*pi*Px())*sin(2*pi*Py())*oneY();
// There gradg is a column vector!
// Use trans() to get a row vector
double normH1_g = normH1( _range=elements(mesh),
_expr=g,
_grad_expr=trans(gradg) );
With test or trial function u
double errorH1 = normH1( _range=elements(mesh),
_expr=(u-g),
_grad_expr=(gradv(u)-trans(gradg)) );
1.3. \(L^\infty\) norm
You can evaluate the infinity norm using the normLinf() function:
1.3.1. Interface
normLinf( _range, _expr, _pset, _geomap );
Required parameters:
-
_range
= domain of integration -
_expr
= mesurable function -
_pset
= set of points (e.g. quadrature points)
Optional parameters:
-
_geomap
= type of geometric mapping.-
Default =
GEOMAP_OPT
-
The normLinf()
function returns not only the maximum of the function
over a sampling of each element thanks to the _pset
argument but
also the coordinates of the point where the function is maximum. The
returned data structure provides the following interface
-
value()
: return the maximum value -
operator()()
: synonym tovalue()
-
arg()
: coordinates of the point where the function is maximum