[dune-pdelab] parallel L2-error norms and parallel volume integrals
Gregor Corbin
corbin at mathematik.uni-kl.de
Tue Dec 11 12:50:24 CET 2018
Hi,
for this i usually use a difference adapter on discrete grid functions.
Probably something similar to this should work:
#include <dune/pdelab/test/l2norm.hh>
#include <dune/pdelab/gridfunctionspace/gridfunctionadapter.hh>
using DGF = Dune::PDELab::DiscreteGridFunction<GFS,RHSType>;
DGF numericDGF(*gfs_p,numeric_sol);
DGF analyticDGF(*gfs_p,analytic_sol);
Dune::PDELab::DifferenceAdapter<DGF,DGF> diff2(numericDGF,analyticDGF);
double l2err2 = l2norm2(diff2);
globall2err = sqrt(gridView.comm().sum(l2err2));
lg,
Gregor Corbin
Am 11.12.18 um 12:03 schrieb Michael Wenske:
> Dear List,
>
> i'm having problems verifying my numerics in a parallel DUNE-PDELab setting.
>
> Im on Dune 2.6, I use a Yasp grid and solve Darcy-like equations.
>
> I constructed an analytical solution, which I can interpolate and write
> to .vtk files
> together with the numerical solution.
> I can compare the numerical solution and the analytic solution with the
> tools
> Paraview provides (Calculator-filter and Integrate Variables etc. ). And
> so far things look wonderful,
> but this way I rely on the correctness of the Paraview implementation.
> And I have to click a lot.
>
> Now I want to do the same thing in the code in a more automated fashion.
> Is there PDELab functionality to -globally- compare and analyze two
> solutions?
>
> I tried:
>
> double get_error_norm(const RHSType& numeric_sol, const RHSType&
> analytic_sol){
> auto difference = numeric_sol-analytic_sol; //Fail
> ISTLHelperType parHelper(*gfs_p);
> Dune::PDELab::OverlappingScalarProduct<GFS,RHSType>
> ovlp_scalar_prod(*gfs_p,parHelper);
> double norm =ovlp_scalar_prod.norm(difference);
> return norm;
> }
>
> where RHSType is :"Dune::PDELab::Backend::Vector<GFS,RangeFieldType>"
>
> but it seems that there are no parallel operators implemented to do *=,
> +, - on these Vector types.
> Technically the above approach is not correct, as not every coefficient
> should be weighted in the same way.
> Coefficients which lie on the Border of the domain (Q1, quadrilateral)
> should be counted half/quarter.
>
> I might force my Local operator to do such operations, since the Methods
> are provided with the necessary
> Information upon assembly, but it seems wrong to me to mix the numerical
> solution steps with the analysis-related
> code.
>
> If I would set up something like a AnalyserLocalOperator, how would I
> traverse this correctly? Currently I pass
> The local operators as template types to the solvers, but in this case
> there is nothing to solve.
> The whole Abstraction of the Local Operator is pointed towards
> assembling a Matrix( for a Solver), so If I just want
> to analyze my solution w.r.t volume, histogram, L2 error against
> analytic sol etc. is this still the correct approach?
>
> This might be trivial for the ol' folks but I can't seem to figure out a
> clean way to do this.
>
> thanks in advance,
>
> Michael Wenske
>
>
>
>
> _______________________________________________
> dune-pdelab mailing list
> dune-pdelab at lists.dune-project.org
> https://lists.dune-project.org/mailman/listinfo/dune-pdelab
-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 819 bytes
Desc: OpenPGP digital signature
URL: <https://lists.dune-project.org/pipermail/dune-pdelab/attachments/20181211/3a54d5fb/attachment.sig>
More information about the dune-pdelab
mailing list