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 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368
|
/*
Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <config.h>
#define BOOST_TEST_MODULE OPM_test_PreconditionerFactory
#include <boost/test/unit_test.hpp>
#include <opm/simulators/linalg/matrixblock.hh>
#include <opm/simulators/linalg/ilufirstelement.hh>
#include <opm/simulators/linalg/PreconditionerFactory_impl.hpp>
#include <opm/simulators/linalg/PropertyTree.hpp>
#include <opm/simulators/linalg/FlexibleSolver.hpp>
#include <opm/simulators/linalg/getQuasiImpesWeights.hpp>
#include <dune/common/fvector.hh>
#include <dune/istl/bvector.hh>
#include <dune/istl/bcrsmatrix.hh>
#include <dune/istl/matrixmarket.hh>
#include <dune/istl/solvers.hh>
#include <fstream>
#include <iostream>
template <class X>
class NothingPreconditioner : public Dune::Preconditioner<X, X>
{
public:
virtual void pre(X&, X&) override
{
}
virtual void apply(X& v, const X& d) override
{
v = d;
}
virtual void post(X&) override
{
}
virtual Dune::SolverCategory::Category category() const override
{
return Dune::SolverCategory::sequential;
}
};
template <int bz>
Dune::BlockVector<Dune::FieldVector<double, bz>>
testPrec(const Opm::PropertyTree& prm, const std::string& matrix_filename, const std::string& rhs_filename)
{
using Matrix = Dune::BCRSMatrix<Opm::MatrixBlock<double, bz, bz>>;
using Vector = Dune::BlockVector<Dune::FieldVector<double, bz>>;
Matrix matrix;
{
std::ifstream mfile(matrix_filename);
if (!mfile) {
throw std::runtime_error("Could not read matrix file");
}
using M = Dune::BCRSMatrix<Dune::FieldMatrix<double, bz, bz>>;
readMatrixMarket(reinterpret_cast<M&>(matrix), mfile); // Hack to avoid hassle
}
Vector rhs;
{
std::ifstream rhsfile(rhs_filename);
if (!rhsfile) {
throw std::runtime_error("Could not read rhs file");
}
readMatrixMarket(rhs, rhsfile);
}
using Operator = Dune::MatrixAdapter<Matrix, Vector, Vector>;
Operator op(matrix);
using PrecFactory = Opm::PreconditionerFactory<Operator,Dune::Amg::SequentialInformation>;
bool transpose = false;
if(prm.get<std::string>("preconditioner.type") == "cprt"){
transpose = true;
}
constexpr int pressure_index = std::min(1, bz-1);
std::function<Vector()> wc = [&matrix, transpose, pressure_index]()
{
return Opm::Amg::getQuasiImpesWeights<Matrix, Vector>(matrix, pressure_index, transpose, false);
};
auto prec = PrecFactory::create(op, prm.get_child("preconditioner"), wc, pressure_index);
Dune::BiCGSTABSolver<Vector> solver(op, *prec, prm.get<double>("tol"), prm.get<int>("maxiter"), prm.get<int>("verbosity"));
Vector x(rhs.size());
Dune::InverseOperatorResult res;
solver.apply(x, rhs, res);
return x;
}
void test1(const Opm::PropertyTree& prm)
{
constexpr int bz = 1;
auto sol = testPrec<bz>(prm, "matr33.txt", "rhs3.txt");
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {-1.62493,
-1.76435e-06,
1.86991e-10,
-458.542,
2.28308e-06,
-2.45341e-07,
-1.48005,
-5.02264e-07,
-1.049e-05};
BOOST_REQUIRE_EQUAL(sol.size(), expected.size());
for (size_t i = 0; i < sol.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(sol[i][row], expected[i][row], 1e-3);
}
}
}
void test3(const Opm::PropertyTree& prm)
{
constexpr int bz = 3;
auto sol = testPrec<bz>(prm, "matr33.txt", "rhs3.txt");
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {{-1.62493, -1.76435e-06, 1.86991e-10},
{-458.542, 2.28308e-06, -2.45341e-07},
{-1.48005, -5.02264e-07, -1.049e-05}};
BOOST_REQUIRE_EQUAL(sol.size(), expected.size());
for (size_t i = 0; i < sol.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(sol[i][row], expected[i][row], 1e-3);
}
}
}
BOOST_AUTO_TEST_CASE(TestDefaultPreconditionerFactory)
{
// Read parameters.
Opm::PropertyTree prm("options_flexiblesolver.json");
// Test with 1x1 block solvers.
test1(prm);
// Test with 3x3 block solvers.
test3(prm);
}
template <int bz>
using M = Dune::BCRSMatrix<Opm::MatrixBlock<double, bz, bz>>;
template <int bz>
using V = Dune::BlockVector<Dune::FieldVector<double, bz>>;
template <int bz>
using O = Dune::MatrixAdapter<M<bz>, V<bz>, V<bz>>;
template <int bz>
using PF = Opm::PreconditionerFactory<O<bz>,Dune::Amg::SequentialInformation>;
BOOST_AUTO_TEST_CASE(TestAddingPreconditioner)
{
// Read parameters.
Opm::PropertyTree prm("options_flexiblesolver_simple.json");
// Test with 1x1 block solvers.
{
constexpr int bz = 1;
BOOST_CHECK_THROW(testPrec<bz>(prm, "matr33.txt", "rhs3.txt"), std::invalid_argument);
}
// Test with 3x3 block solvers.
{
constexpr int bz = 3;
BOOST_CHECK_THROW(testPrec<bz>(prm, "matr33.txt", "rhs3.txt"), std::invalid_argument);
}
// Add preconditioner to factory for block size 1.
PF<1>::addCreator("nothing", [](const O<1>&, const Opm::PropertyTree&, const std::function<V<1>()>&,
std::size_t) {
return Dune::getDummyUpdateWrapper<NothingPreconditioner<V<1>>>();
});
// Test with 1x1 block solvers.
test1(prm);
// Test with 3x3 block solvers.
{
constexpr int bz = 3;
BOOST_CHECK_THROW(testPrec<bz>(prm, "matr33.txt", "rhs3.txt"), std::invalid_argument);
}
// Add preconditioner to factory for block size 3.
PF<3>::addCreator("nothing", [](const O<3>&, const Opm::PropertyTree&, const std::function<V<3>()>&,
std::size_t) {
return Dune::getDummyUpdateWrapper<NothingPreconditioner<V<3>>>();
});
// Test with 1x1 block solvers.
test1(prm);
// Test with 3x3 block solvers.
test3(prm);
}
template<class Mat, class Vec>
class RepeatingOperator : public Dune::AssembledLinearOperator<Mat, Vec, Vec>
{
public:
using matrix_type = Mat;
using domain_type = Vec;
using range_type = Vec;
using field_type = typename Vec::field_type;
Dune::SolverCategory::Category category() const override
{
return Dune::SolverCategory::sequential;
}
RepeatingOperator(const Mat& matrix, const int repeats)
: matrix_(matrix)
, repeats_(repeats)
{
}
// y = A*x;
virtual void apply(const Vec& x, Vec& y) const override
{
y = 0;
applyscaleadd(1.0, x, y);
}
// y += \alpha * A * x
virtual void applyscaleadd(field_type alpha, const Vec& x, Vec& y) const override
{
Vec temp1 = x;
Vec temp2 = x; // For size.
temp2 = 0.0;
for (int rr = 0; rr < repeats_; ++rr) {
// mv below means: temp2 = matrix_ * temp1;
matrix_.mv(temp1, temp2);
temp1 = temp2;
}
temp2 *= alpha;
y += temp2;
}
virtual const matrix_type& getmat() const override
{
return matrix_;
}
protected:
const Mat& matrix_;
const int repeats_;
};
template <int bz>
Dune::BlockVector<Dune::FieldVector<double, bz>>
testPrecRepeating(const Opm::PropertyTree& prm, const std::string& matrix_filename, const std::string& rhs_filename)
{
using Matrix = M<bz>;
using Vector = V<bz>;
Matrix matrix;
{
std::ifstream mfile(matrix_filename);
if (!mfile) {
throw std::runtime_error("Could not read matrix file");
}
using M = Dune::BCRSMatrix<Dune::FieldMatrix<double, bz, bz>>;
#pragma GCC diagnostic push
#pragma GCC diagnostic ignored "-Wstrict-aliasing"
readMatrixMarket(reinterpret_cast<M&>(matrix), mfile); // Hack to avoid hassle
#pragma GCC diagnostic pop
}
Vector rhs;
{
std::ifstream rhsfile(rhs_filename);
if (!rhsfile) {
throw std::runtime_error("Could not read rhs file");
}
readMatrixMarket(rhs, rhsfile);
}
using Operator = RepeatingOperator<Matrix, Vector>;
Operator op(matrix, 2);
using PrecFactory = Opm::PreconditionerFactory<Operator,Dune::Amg::SequentialInformation>;
// Add no-oppreconditioner to factory for block size 1.
PrecFactory::addCreator("nothing", [](const Operator&, const Opm::PropertyTree&, const std::function<Vector()>&,
std::size_t) {
return Dune::getDummyUpdateWrapper<NothingPreconditioner<Vector>>();
});
auto prec = PrecFactory::create(op, prm.get_child("preconditioner"));
Dune::BiCGSTABSolver<Vector> solver(op, *prec, prm.get<double>("tol"), prm.get<int>("maxiter"), prm.get<int>("verbosity"));
Vector x(rhs.size());
Dune::InverseOperatorResult res;
solver.apply(x, rhs, res);
return x;
}
void test1rep(const Opm::PropertyTree& prm)
{
constexpr int bz = 1;
auto sol = testPrecRepeating<bz>(prm, "matr33rep.txt", "rhs3rep.txt");
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {0.285714285714286,
0.285714285714286,
0.285714285714286,
-0.214285714285714,
-0.214285714285714,
-0.214285714285714,
-0.214285714285714,
-0.214285714285714,
-0.214285714285714};
BOOST_REQUIRE_EQUAL(sol.size(), expected.size());
for (size_t i = 0; i < sol.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(sol[i][row], expected[i][row], 1e-3);
}
}
}
void test3rep(const Opm::PropertyTree& prm)
{
constexpr int bz = 3;
auto sol = testPrecRepeating<bz>(prm, "matr33rep.txt", "rhs3rep.txt");
Dune::BlockVector<Dune::FieldVector<double, bz>> expected {
{0.285714285714286, 0.285714285714286, 0.285714285714286},
{-0.214285714285714, -0.214285714285714, -0.214285714285714},
{-0.214285714285714, -0.214285714285714, -0.214285714285714}
};
BOOST_REQUIRE_EQUAL(sol.size(), expected.size());
for (size_t i = 0; i < sol.size(); ++i) {
for (int row = 0; row < bz; ++row) {
BOOST_CHECK_CLOSE(sol[i][row], expected[i][row], 1e-3);
}
}
}
BOOST_AUTO_TEST_CASE(TestWithRepeatingOperator)
{
// Read parameters.
Opm::PropertyTree prm("options_flexiblesolver_simple.json");
// Test with 1x1 block solvers.
test1rep(prm);
// Test with 3x3 block solvers.
test3rep(prm);
}
|