DUNE PDELab (git)
Operations on grid functions
A frequent example of pointwise operations on GridFunctions is calculating error norms (e.g. analytical vs numerical solution). Here we take the squared difference of two functions. However, the same approach applies to any operation.
First, let's define some grid functions for testing.
auto analyticFunction1 = Dune::PDELab::makeGridFunctionFromCallable (grid.leafGridView(), [&](const auto& x){
return exp(-(x*x));
});
auto analyticFunction2 = Dune::PDELab::makeGridFunctionFromCallable (grid.leafGridView(), [&](const auto& x){
return x[0]*3.0;
});
Now, we can either plug them into an existing adapter (itself a GridFunction) computing the difference squared
Dune::PDELab::DifferenceSquaredAdapter<decltype(analyticFunction1),decltype(analyticFunction2)> difference(analyticFunction1, analyticFunction2);
// Note: Newer compiler versions do not require template arguments here!
Adapter returning ||f1(x)-f2(x)||^2 for two given grid functions.
Definition: gridfunctionadapter.hh:70
or define it ourselves. Here we use the per-element evaluation interface receiving an element and local coordinates on the reference element.
auto difference = Dune::PDELab::makeGridFunctionFromCallable (grid.leafGridView(), [&](const auto& element, const auto& xlocal){
Dune::FieldVector<double,1> y1;
analyticFunction1.evaluate(element, xlocal, y1);
Dune::FieldVector<double,1> y2;
analyticFunction2.evaluate(element, xlocal, y2);
y2 -= y1;
return y2.two_norm2();
});
We can use the resulting GridFunctions as always, e.g. integrate them:
auto integral = Dune::PDELab::integrateGridFunction(difference,10);
GF::Traits::RangeType integrateGridFunction(const GF &gf, unsigned qorder=1)
Integrate a GridFunction.
Definition: functionutilities.hh:51
Full example code: recipe-grid-function-operations.cc
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)