convergence.cc 7.83 KB
Newer Older
1
2
3
4
5
6
7
8
9
// -*- 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 <dune/common/parallel/mpihelper.hh> // An initializer of MPI
#include <dune/curvedsurfacegrid/curvedsurfacegrid.hh>
10
#include <dune/curvedsurfacegrid/gridfunction/analyticgridfunction.hh>
11
12
13
14
15
16
#include <dune/geometry/quadraturerules.hh>
#include <dune/grid/io/file/gmshreader.hh>
#include <dune/localfunctions/lagrange/pk.hh>
//#include <dune/grid/albertagrid.hh>
#include <dune/foamgrid/foamgrid.hh>

17
#include <dune/functions/common/differentiablefunctionfromcallables.hh>
18
//#include <dune/functions/gridfunctions/analyticgridviewfunction.hh>
19

20
21
22
#define STR(s) STR_HELPER(s)
#define STR_HELPER(s) #s

23
const int order = 4;
24
const int quad_order = order+10;
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
const int num_levels = 5;

// Test Hausdorf-distance
template <class Grid, class Projection>
typename Grid::ctype inf_error(Grid const& grid, Projection const& projection)
{
  using QuadProvider = Dune::QuadratureRules<typename Grid::ctype, 2>;
  typename Grid::ctype dist = 0;
  for (auto const& element : elements(grid.leafGridView()))
  {
    auto geometry = element.geometry();

    auto const& quadRule = QuadProvider::rule(element.type(), quad_order);
    for (auto const& qp : quadRule) {
      auto x = geometry.global(qp.position());
      auto y = projection(x);
      dist = std::max(dist, (x - y).two_norm());
    }
  }

  return dist;
}

// Test integrated Hausdorf-distance
template <class Grid, class Projection>
typename Grid::ctype L2_error(Grid const& grid, Projection const& projection)
{
  using QuadProvider = Dune::QuadratureRules<typename Grid::ctype, 2>;
  typename Grid::ctype dist = 0;
  for (auto const& element : elements(grid.leafGridView()))
  {
    auto geometry = element.geometry();

    auto const& quadRule = QuadProvider::rule(element.type(), quad_order);
    for (auto const& qp : quadRule) {
      auto x = geometry.global(qp.position());
      auto y = projection(x);
      dist += (x - y).two_norm2() * qp.weight() * geometry.integrationElement(qp.position());
    }
  }

  return std::sqrt(dist);
}

// Test integrated error in normal vectors
template <class Grid, class Projection>
typename Grid::ctype normal_error(Grid const& grid, Projection const& projection)
{
  using QuadProvider = Dune::QuadratureRules<typename Grid::ctype, 2>;
  typename Grid::ctype dist = 0;
  for (auto const& element : elements(grid.leafGridView()))
  {
    auto geometry = element.geometry();
    auto const& quadRule = QuadProvider::rule(element.type(), quad_order);
    for (auto const& qp : quadRule) {
      auto x = geometry.impl().normal(qp.position());
      auto y = projection(geometry.global(qp.position()));
      if ((y + x).two_norm() < 0.9)
        x *= -1;

      dist += (x - y).two_norm2() * qp.weight() * geometry.integrationElement(qp.position());
    }
  }

  return std::sqrt(dist);
}

// Test integrated error in mean curvature
template <class Grid, class Projection>
typename Grid::ctype curvature_error(Grid const& grid, Projection const& projection)
{
  using GlobalCoordinate = typename Grid::template Codim<0>::Entity::Geometry::GlobalCoordinate;
  using QuadProvider = Dune::QuadratureRules<typename Grid::ctype, 2>;
98
  using FiniteElementType = Dune::PkLocalFiniteElement<typename Grid::ctype, typename Grid::ctype, 2, order+4>;
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
  FiniteElementType fe;

  using LBTraits = typename FiniteElementType::Traits::LocalBasisType::Traits;
  std::vector<typename LBTraits::JacobianType> shapeGradients;
  std::vector<GlobalCoordinate> gradients;
  std::vector<GlobalCoordinate> normals;

  typename Grid::ctype dist = 0;
  for (auto const& element : elements(grid.leafGridView()))
  {
    auto geometry = element.geometry();
    auto const& localBasis = fe.localBasis();
    auto const& localInterpolation = fe.localInterpolation();

    auto f = [&](auto const& local) -> GlobalCoordinate
    {
      auto x = geometry.impl().normal(local);
      auto y = projection(geometry.global(local));
      if ((y + x).two_norm() < 0.9)
        x *= -1;
      return x;
    };
    localInterpolation.interpolate(f, normals);

    shapeGradients.resize(localBasis.size());
    gradients.resize(localBasis.size());

    auto const& quadRule = QuadProvider::rule(element.type(), quad_order);
    for (auto const& qp : quadRule) {
      auto jInvT = geometry.jacobianInverseTransposed(qp.position());
      localBasis.evaluateJacobian(qp.position(), shapeGradients);

      for (std::size_t i = 0; i < gradients.size(); ++i)
        jInvT.mv(shapeGradients[i][0], gradients[i]);

      typename Grid::ctype H = 0;
      for (std::size_t i = 0; i < gradients.size(); ++i)
        H += normals[i].dot(gradients[i]);

      dist += std::abs(H - 2) * qp.weight() * geometry.integrationElement(qp.position());
    }
  }

  return dist;
}

// longest edge in the grid
template <class Grid>
typename Grid::ctype edge_length(Grid const& grid)
{
  typename Grid::ctype h = 0;
150
  for (auto const& e : edges(grid.hostGrid().leafGridView()))
151
152
153
154
155
    h = std::max(h, e.geometry().volume());

  return h;
}

156

157
158
159
int main(int argc, char** argv)
{
  using namespace Dune;
160
  MPIHelper::instance(argc, argv);
161
162
163
164
165

  // using HostGrid = AlbertaGrid<2,3>;
  using HostGrid = FoamGrid<2,3>;
  std::unique_ptr<HostGrid> hostGrid = GmshReader<HostGrid>::read( STR(DUNE_GRID_PATH) "sphere_coarse.msh");

166
167
168
169
170
171
172
  using Signature = FieldVector<double,3>(FieldVector<double,3>);
  auto sphere = Functions::makeDifferentiableFunctionFromCallables(
    Functions::SignatureTag<Signature>{},
    [](auto const& x) {
      return x / x.two_norm();
    },
    [](auto const& x) {
173
      double nrm = x.two_norm();
174
175
176
      FieldMatrix<double,3,3> out;
      for (int i = 0; i < 3; ++i)
        for (int j = 0; j < 3; ++j)
177
          out[i][j] = ((i == j ? 1 : 0) - (x[i]/nrm) * (x[j]/nrm)) / nrm;
178
179
180
      return out;
    }
  );
181
182
  //auto sphereGridFct = Functions::makeAnalyticGridViewFunction(sphere, hostGrid->leafGridView());
  auto sphereGridFct = analyticGridFunction<HostGrid>(sphere);
183
184
185

  using Grid = CurvedSurfaceGrid<decltype(sphereGridFct)>;
  Grid grid(*hostGrid, sphereGridFct);
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210

  std::vector<typename Grid::ctype> inf_errors, L2_errors, normal_errors, curvature_errors, edge_lengths;
  for (int i = 0; i < num_levels; ++i) {
    if (i > 0)
      grid.globalRefine(1);

    inf_errors.push_back( inf_error(grid, sphere) );
    L2_errors.push_back( L2_error(grid, sphere) );
    normal_errors.push_back( normal_error(grid, sphere) );
    curvature_errors.push_back( curvature_error(grid, sphere) );
    edge_lengths.push_back( edge_length(grid) );
  }


  std::vector<typename Grid::ctype> eocInf(num_levels, 0), eocL2(num_levels, 0), eocNormal(num_levels, 0), eocCurvature(num_levels, 0);
  for (int i = 1; i < num_levels; ++i) {
    eocInf[i]       = std::log(inf_errors[i]/inf_errors[i-1]) / std::log(edge_lengths[i]/edge_lengths[i-1]);
    eocL2[i]        = std::log(L2_errors[i]/L2_errors[i-1]) / std::log(edge_lengths[i]/edge_lengths[i-1]);
    eocNormal[i]    = std::log(normal_errors[i]/normal_errors[i-1]) / std::log(edge_lengths[i]/edge_lengths[i-1]);
    eocCurvature[i] = std::log(curvature_errors[i]/curvature_errors[i-1]) / std::log(edge_lengths[i]/edge_lengths[i-1]);
  }

  auto print_line = [](auto i, auto const& data) {
    std::cout << i;
    for (auto const& d : data)
211
      std::cout << '\t' << "| " << d;
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
    std::cout << std::endl;
  };

  auto print_break = [] {
    std::cout.width(8 + 8*16);
    std::cout.fill('-');
    std::cout << '-' << std::endl;
    std::cout.fill(' ');
  };

  std::cout.setf(std::ios::scientific);
  print_line("level", std::vector<std::string>{"err_inf", "eoc_inf", "err_L2", "eoc_L2", "err_norm", "eoc_norm", "err_curv", "eoc_curv"});
  print_break();

  for (int i = 0; i < num_levels; ++i)
    print_line(i, std::vector<double>{inf_errors[i], eocInf[i], L2_errors[i], eocL2[i], normal_errors[i], eocNormal[i], curvature_errors[i], eocCurvature[i]});
}