[dune-functions] Problems writing vector-valued vtk files
Steffen Müthing
steffen.muething at iwr.uni-heidelberg.de
Wed Mar 18 00:06:11 CET 2015
> Am 17.03.2015 um 15:32 schrieb Oliver Sander <sander at igpm.rwth-aachen.de>:
>
> Hi guys,
> I am trying to get the new vtkwriter to write vector-valued functions,
> and I can't get it to work. Is this a bug or am I forgetting some important
> option? To reproduce, please apply the attached short patch to poisson-pq2.cc.
> Then, running the program aborts with
>
> poisson-pq2: /home/sander/dune/dune-common/dune/common/densevector.hh:538: typename Dune::PromotionTraits<typename Dune::DenseMatVecTraits<V>::value_type, typename
> Dune::DenseVector<Other>::field_type>::PromotedType Dune::DenseVector<V>::operator*(const Dune::DenseVector<Other>&) const [with Other = Dune::FieldVector<double, 1>; V = Dune::FieldVector<double, 3>;
> typename Dune::PromotionTraits<typename Dune::DenseMatVecTraits<V>::value_type, typename Dune::DenseVector<Other>::field_type>::PromotedType = double]: Assertion `y.size() == size()' failed.
> Abgebrochen
That’s strange - I can run this program and get a valid VTK file. Are you sure your dune-grid
is up to date?
Steffen
>
> Thanks for your help,
> Oliver
>
>
>
> diff --git a/examples/poisson-pq2.cc b/examples/poisson-pq2.cc
> index 944aacc..6b5a079 100644
> --- a/examples/poisson-pq2.cc
> +++ b/examples/poisson-pq2.cc
> @@ -20,7 +20,7 @@
> #include <dune/istl/solvers.hh>
>
> #include <dune/functions/functionspacebases/interpolate.hh>
> -#include <dune/functions/functionspacebases/pq2nodalbasis.hh>
> +#include <dune/functions/functionspacebases/pq1nodalbasis.hh>
> #include <dune/functions/gridfunctions/discretescalarglobalbasisfunction.hh>
> #include <dune/functions/gridfunctions/gridviewfunction.hh>
>
> @@ -304,7 +304,7 @@ int main (int argc, char *argv[]) try
> // Choose a finite element space
> /////////////////////////////////////////////////////////
>
> - typedef Functions::PQ2NodalBasis<GridView> FEBasis;
> + typedef Functions::PQ1NodalBasis<GridView> FEBasis;
> FEBasis feBasis(gridView);
>
> /////////////////////////////////////////////////////////
> @@ -399,6 +399,7 @@ int main (int argc, char *argv[]) try
> // Make a discrete function from the FE basis and the coefficient vector
> ////////////////////////////////////////////////////////////////////////////
>
> +#if 0
> Dune::Functions::DiscreteScalarGlobalBasisFunction<decltype(feBasis),decltype(x)> xFunction(feBasis,x);
> auto localXFunction = localFunction(xFunction);
>
> @@ -410,6 +411,23 @@ int main (int argc, char *argv[]) try
> SubsamplingVTKWriter<GridView> vtkWriter(gridView,2);
> vtkWriter.addVertexData(localXFunction, VTK::FieldInfo("x", VTK::FieldInfo::Type::scalar, 1));
> vtkWriter.write("poisson-pq2");
> +#else
> + typedef BlockVector<FieldVector<double,3> > EmbeddedVectorType;
> + EmbeddedVectorType xEmbedded(x.size());
> + for (size_t i=0; i<x.size(); i++) {
> + xEmbedded[i][0] = 0;
> + xEmbedded[i][1] = 0;
> + xEmbedded[i][2] = 0;
> + }
> +
> + Dune::Functions::DiscreteScalarGlobalBasisFunction<decltype(feBasis),decltype(xEmbedded)> xFunction(feBasis,xEmbedded);
> + auto localXFunction = localFunction(xFunction);
> +
> + VTKWriter<GridView> vtkWriter(gridView);
> + vtkWriter.addVertexData(localXFunction, VTK::FieldInfo("x", VTK::FieldInfo::Type::scalar, 3));
> +
> + vtkWriter.write("test");
> +#endif
>
> }
> // Error handling
> <provoke_vector_vtk_error.patch>_______________________________________________
> dune-functions mailing list
> dune-functions at dune-project.org
> http://lists.dune-project.org/mailman/listinfo/dune-functions
-------------- next part --------------
A non-text attachment was scrubbed...
Name: signature.asc
Type: application/pgp-signature
Size: 495 bytes
Desc: Message signed with OpenPGP using GPGMail
URL: <https://lists.dune-project.org/pipermail/dune-functions/attachments/20150317/843b1f3c/attachment.sig>
More information about the dune-functions
mailing list