Chucuim
> If I don't have the exact solution, but I want to compute the L2 form of
> u_{h}-u_{finest} (i.e. see u_{finest} as the exact solution, and
> u_{finest} is the numerical solution with the finest mesh grid among the
> series we want to test, for example, the width of the unit meshgrid is h
> = 1/2^k, k=6,7,8,9,10. So the u_{finest}=u_{h=1/2^10}). Does deal.ii
> have a function like
> |
> template<intdim,classInVector,classOutVector,intspacedim>
>
> voidVectorTools::integrate_difference ( constMapping
> <https://www.dealii.org/developer/doxygen/deal.II/classMapping.html><dim,spacedim
>
> >& mapping,
> constDoFHandler
> <https://www.dealii.org/developer/doxygen/deal.II/classDoFHandler.html><dim,spacedim
>
> >& dof,
> constInVector& fe_function,
> constFunction
> <https://www.dealii.org/developer/doxygen/deal.II/classFunction.html><spacedim,double>&
>
> exact_solution,
> OutVector& difference,
> constQuadrature
> <https://www.dealii.org/developer/doxygen/deal.II/classQuadrature.html><dim
> >& q,
> constNormType
> <https://www.dealii.org/developer/doxygen/deal.II/group__numerics.html#ga69967cb7a148a7169963126249213db1>&
>
> norm,
> constFunction
> <https://www.dealii.org/developer/doxygen/deal.II/classFunction.html><spacedim,double>*
>
> weight =|nullptr|,
> constdouble exponent =|2.|
> )
>
> |
>
>
> and doesn't use
>
> constFunction
> <https://www.dealii.org/developer/doxygen/deal.II/classFunction.html><
> spacedim,double>&
> exact_solution,
>
>
> but use
>
> constInVector& fe_function,
>
> at this location of input parameters. Thanks in advance!
You can't do this directly -- the function just wants an object that can
be evaluated at arbitrary points. Passing a vector by itself will not
help: a vector is just a bunch of numbers, but it does not by itself
explain how it is to be interepreted. Indeed, a vector by itself is
meaningless unless you also specify a DoFHandler that makes the
connection between nodal values, where these nodes are located in
physical space, and how the shape functions are defined.
But there is a class that combines all of this knowledge and that then
represents a finite element field (=DoFHandler + Vector) as a function:
FEFieldFunction. Take a look at that class, as it implements the
interface of a Function object that can be passed in the place you point
out. It is not particularly fast, but should allow you to do what you want.
Of course, you have to make sure that both the 'dof' and 'fe_function'
arguments to the VectorTools::integrate_difference() function match
(i.e., you can't just save a solution vector, and then refine the
triangulation and thereby change the DoFHandler), as well as the
DoFHandler and Vector objects you pass to FEFieldFunction.
Best
W.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
The deal.II project is located at http://www.dealii.org/
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
---
You received this message because you are subscribed to the Google Groups
"deal.II User Group" group.
To unsubscribe from this group and stop receiving emails from it, send an email
to [email protected].
For more options, visit https://groups.google.com/d/optout.