vtkwriter.cc 2.78 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:

#ifdef HAVE_CONFIG_H
# include "config.h"
#endif

#include <iostream>
#include <vector>

#include <dune/common/parallel/mpihelper.hh> // An initializer of MPI
#include <dune/common/exceptions.hh> // We use exceptions
#include <dune/common/filledarray.hh>

#include <dune/functions/functionspacebases/defaultglobalbasis.hh>
#include <dune/functions/functionspacebases/lagrangebasis.hh>
#include <dune/functions/functionspacebases/interpolate.hh>
#include <dune/functions/gridfunctions/analyticgridviewfunction.hh>
#include <dune/functions/gridfunctions/discreteglobalbasisfunction.hh>
#include <dune/grid/uggrid.hh>
#include <dune/grid/yaspgrid.hh>

23
#include <dune/vtk/writers/vtkunstructuredgridwriter.hh>
Praetorius, Simon's avatar
Praetorius, Simon committed
24
25

using namespace Dune;
26
using namespace Dune::experimental;
Praetorius, Simon's avatar
Praetorius, Simon committed
27
28
using namespace Dune::Functions;

Praetorius, Simon's avatar
Praetorius, Simon committed
29
#define GRID_TYPE 1
Praetorius, Simon's avatar
Praetorius, Simon committed
30
31
32
33
34
35
36
37
38
39

int main(int argc, char** argv)
{
  Dune::MPIHelper::instance(argc, argv);

  const int dim = 3;

#if GRID_TYPE == 1
  using GridType = YaspGrid<dim>;
  FieldVector<double,dim> upperRight; upperRight = 1.0;
Praetorius, Simon's avatar
Praetorius, Simon committed
40
  auto numElements = filledArray<dim,int>(8);
Praetorius, Simon's avatar
Praetorius, Simon committed
41
42
43
44
45
  GridType grid(upperRight,numElements);
#elif GRID_TYPE == 2
  using GridType = UGGrid<dim>;
  FieldVector<double,dim> lowerLeft; lowerLeft = 0.0;
  FieldVector<double,dim> upperRight; upperRight = 1.0;
46
  auto numElements = filledArray<dim,unsigned int>(4);
Praetorius, Simon's avatar
Praetorius, Simon committed
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
  auto gridPtr = StructuredGridFactory<GridType>::createSimplexGrid(lowerLeft, upperRight, numElements);
  auto& grid = *gridPtr;
#endif

  using GridView = typename GridType::LeafGridView;
  GridView gridView = grid.leafGridView();

  using namespace BasisFactory;
  auto basis = makeBasis(gridView, lagrange<1>());

  std::vector<double> p1function(basis.dimension());

  interpolate(basis, p1function, [](auto const& x) {
    return 100*x[0] + 10*x[1] + 1*x[2];
  });

  // write discrete global-basis function
  auto p1FctWrapped = makeDiscreteGlobalBasisFunction<double>(basis, p1function);

Praetorius, Simon's avatar
Praetorius, Simon committed
66
  using Writer = VtkUnstructuredGridWriter<GridView>;
Praetorius, Simon's avatar
Praetorius, Simon committed
67
  Writer vtkWriter(gridView);
68
  vtkWriter.addPointData(p1FctWrapped, "p1");
Praetorius, Simon's avatar
Praetorius, Simon committed
69
70
71
72
73
74
75
  vtkWriter.addCellData(p1FctWrapped, "p0");

  // write analytic function
  auto p1Analytic = makeAnalyticGridViewFunction([](auto const& x) {
    return std::sin(10*x[0])*std::cos(10*x[1])+std::sin(10*x[2]);
  }, gridView);

76
  vtkWriter.addPointData(p1Analytic, "analytic");
Praetorius, Simon's avatar
Praetorius, Simon committed
77
78
79
80
81
82
83
84

  vtkWriter.write("test_ascii_float32.vtu", Vtk::ASCII);
  vtkWriter.write("test_binary_float32.vtu", Vtk::BINARY);
  vtkWriter.write("test_compressed_float32.vtu", Vtk::COMPRESSED);
  vtkWriter.write("test_ascii_float64.vtu", Vtk::ASCII, Vtk::FLOAT64);
  vtkWriter.write("test_binary_float64.vtu", Vtk::BINARY, Vtk::FLOAT64);
  vtkWriter.write("test_compressed_float64.vtu", Vtk::COMPRESSED, Vtk::FLOAT64);
}