File: istlvectorbackendtest.cc

package info (click to toggle)
dune-functions 2.10.0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 1,544 kB
  • sloc: cpp: 14,241; python: 661; makefile: 3
file content (220 lines) | stat: -rw-r--r-- 7,118 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
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
98
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
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
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
211
212
213
214
215
216
217
218
219
220
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:

// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later

#include <config.h>

#include <vector>
#include <array>

#include <dune/common/parallel/mpihelper.hh>
#include <dune/common/reservedvector.hh>
#include <dune/common/classname.hh>
#include <dune/common/fvector.hh>
#include <dune/common/indices.hh>
#include <dune/common/tuplevector.hh>
#include <dune/common/test/testsuite.hh>

#include <dune/istl/bvector.hh>
#include <dune/istl/multitypeblockvector.hh>

#include <dune/functions/common/type_traits.hh>
#include <dune/functions/backends/concepts.hh>
#include <dune/functions/backends/istlvectorbackend.hh>


using namespace Dune;



/**
 * \brief A Dummy global basis
 *
 * This is a mock class providing non-uniform size information.
 * It's non-uniform in the sense that not all multi-indices
 * have the same length.
 */
template<std::size_t dim>
class GlobalBasisMoc
{
public:
    using size_type = std::size_t;
    using SizePrefix = Dune::ReservedVector<std::size_t, 3>;
    using MultiIndex = Dune::ReservedVector<std::size_t, 3>;

    /**
     * \brief Construct from basis
     */
    GlobalBasisMoc()
    {}

    /**
     * \brief Return number possible values for next position in multi index
     *
     * This shall vanish. It's just here such that this can be used
     * as size provider n place of the basis.
     */
    size_type size(const SizePrefix& prefix) const
    {
      if (prefix.size() == 0)
        return 2;
      if (prefix.size() == 1)
      {
        if (prefix[0] == 0)
          return 23;
        if (prefix[0] == 1)
          return 42;
      }
      if (prefix.size() == 2)
      {
        if (prefix[0] == 0)
          return dim;
        if (prefix[0] == 1)
          return 0;
      }
      if (prefix.size() == 3)
        return 0;
      return 0;
    }

    operator size_type () const
    {
        return 23*dim+42;
    }

};


template<class Vector, class Coefficient, std::size_t dim, class MultiIndex>
Dune::TestSuite checkISTLVectorBackend(std::string shortName="")
{
  Dune::TestSuite test(shortName);

  using namespace Dune::Indices;
  using Basis = GlobalBasisMoc<dim>;
  using SizePrefix = typename Basis::SizePrefix;

  Basis basis;

  // Create raw vector
  Vector x_raw;

  // Create wrapped vector
  auto x = Dune::Functions::istlVectorBackend(x_raw);

  test.require(Dune::models<Dune::Functions::Concept::VectorBackend<Basis>, decltype(x)>(), "VectorBackend concept check")
    << "Object returned by istlVectorBackend() does not model the VectorBackend concept";

  // Resize wrapped vector using basis
  x.resize(basis);

  // Derive size information from vector
  test.require(x_raw.size() == basis.size(SizePrefix{}), "resize check")
    << "x_raw.size() is " << x_raw.size() << " but should be " << basis.size(SizePrefix{});

  test.require(x_raw[_0].size() == basis.size(SizePrefix{0}), "resize check")
    << "x_raw[_0].size() is " << x_raw[_0].size() << " but should be " << basis.size(SizePrefix{0});

  for (std::size_t i=0; i<basis.size({0}); ++i)
    test.require(x_raw[_0][i].size() == basis.size(SizePrefix{0,i}), "resize check")
      << "x_raw[_0][" << i << "].size() is " << x_raw[_0][i].size() << " but should be " << basis.size(SizePrefix{0,i});

  test.require(x_raw[_1].size() == basis.size(SizePrefix{1}), "resize check")
    << "x_raw[_1].size() is " << x_raw[_0].size() << " but should be " << basis.size(SizePrefix{1});


  // Assign values to each vector entry
  for (std::size_t i=0; i<x_raw[_0].size(); ++i)
    for (std::size_t j=0; j<x_raw[_0][i].size(); ++j)
      x[MultiIndex{{0, i, j}}] = 0+i+j;
  for (std::size_t i=0; i<x_raw[_1].size(); ++i)
    x[MultiIndex{{1, i}}] = 1+i;


  // Access vector entries via const reference
  const auto& x_const = x;
  for (std::size_t i=0; i<x_raw[_0].size(); ++i)
    for (std::size_t j=0; j<x_raw[_0][i].size(); ++j)
    {
      test.check(x_const[MultiIndex{{0, i, j}}] == Coefficient(0+i+j))
        << "x[{0," << i << "," << j << "}] contains wrong value";
    }
  for (std::size_t i=0; i<x_raw[_1].size(); ++i)
  {
    test.check(x_const[MultiIndex{{1, i}}] == Coefficient(1+i))
      << "x[{1," << i << "}] contains wrong value";
  }

  return test;
}


int main (int argc, char *argv[]) try
{
  // Set up MPI, if available
  MPIHelper::instance(argc, argv);

  ///////////////////////////////////
  //   Generate the grid
  ///////////////////////////////////


  Dune::TestSuite test;
  {
    using VelocityVector = std::vector<std::vector<double>>;
    using PressureVector = std::vector<double>;
    using Coefficient = double;
    using Vector = Dune::TupleVector<VelocityVector, PressureVector>;
    using MultiIndex = ReservedVector<std::size_t, 3>;
    test.subTest(checkISTLVectorBackend<Vector, Coefficient, 2, MultiIndex>("TV<V<V<double>>, V<double>>"));
  }

  {
    using VelocityVector = std::vector<Dune::BlockVector<Dune::FieldVector<double,1>>>;
    using PressureVector = std::vector<Dune::FieldVector<double,1>>;
    using Coefficient = double;
    using Vector = Dune::TupleVector<VelocityVector, PressureVector>;
    using MultiIndex = ReservedVector<std::size_t, 3>;
    test.subTest(checkISTLVectorBackend<Vector, Coefficient, 2, MultiIndex>("TV<V<BV<FV<double,1>>>, V<FV<doule,1>>>"));
  }

  {
    static const std::size_t dim = 5;
    using VelocityVector = std::vector<std::array<Dune::FieldVector<double,1>,dim>>;
    using PressureVector = std::vector<double>;
    using Coefficient = double;
    using Vector = Dune::TupleVector<VelocityVector, PressureVector>;
    using MultiIndex = ReservedVector<std::size_t, 3>;
    test.subTest(checkISTLVectorBackend<Vector, Coefficient, dim, MultiIndex>("TV<V<A<FV<double,1>,5>>, V<double>>"));
  }

  {
    static const std::size_t dim = 5;
    using VelocityVector = Dune::BlockVector<Dune::FieldVector<double,dim>>;
    using PressureVector = Dune::BlockVector<Dune::FieldVector<double,1>>;
    using Coefficient = double;
    using Vector = Dune::MultiTypeBlockVector<VelocityVector, PressureVector>;
    using MultiIndex = ReservedVector<std::size_t, 3>;
    test.subTest(checkISTLVectorBackend<Vector, Coefficient, dim, MultiIndex>("MTBV<BV<FV<double,5>>, BV<FV<double,1>>>"));
  }

  {
    static const std::size_t dim = 3;
    using VelocityVector = std::vector<Dune::MultiTypeBlockVector<Dune::FieldVector<double,1>, double, Dune::FieldVector<double,1>>>;
    using PressureVector = Dune::BlockVector<Dune::FieldVector<double,1>>;
    using Coefficient = double;
    using Vector = Dune::MultiTypeBlockVector<VelocityVector, PressureVector>;
    using MultiIndex = ReservedVector<std::size_t, 3>;
    test.subTest(checkISTLVectorBackend<Vector, Coefficient, dim, MultiIndex>("MTBV<V<MTBV<FV<double,1>, double, FV<double,1>>>, BV<FV<double,1>>"));
  }


  return test.exit();
}
// Error handling
catch (Exception& e) {
  std::cout << e.what() << std::endl;
  return 1;
}