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
|
/*
Copyright 2018 SINTEF Digital, Mathematics and Cybernetics.
Copyright 2018 Equinor.
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 TestRstConv
#define BOOST_TEST_NO_MAIN
#include <boost/test/unit_test.hpp>
#include <boost/version.hpp>
#if BOOST_VERSION / 100000 == 1 && BOOST_VERSION / 100 % 1000 > 66
#include <boost/test/data/test_case.hpp>
#endif
#include <dune/common/parallel/mpihelper.hh>
#include <dune/common/fvector.hh>
#include <dune/istl/bvector.hh>
#include <opm/input/eclipse/Schedule/RSTConfig.hpp>
#include <opm/simulators/flow/RSTConv.hpp>
#include <algorithm>
#include <cstddef>
#include <numeric>
#include <random>
#if HAVE_MPI
struct MPIError
{
MPIError(std::string s, int e) : errorstring(std::move(s)), errorcode(e){}
std::string errorstring;
int errorcode;
};
void MPI_err_handler(MPI_Comm*, int* err_code, ...)
{
std::vector<char> err_string(MPI_MAX_ERROR_STRING);
int err_length;
MPI_Error_string(*err_code, err_string.data(), &err_length);
std::string s(err_string.data(), err_length);
std::cerr << "An MPI Error ocurred:" << std::endl << s << std::endl;
throw MPIError(s, *err_code);
}
#endif
bool
init_unit_test_func()
{
return true;
}
struct TestCase {
std::size_t N;
std::array<int,6> phase;
};
std::ostream& operator<<(std::ostream& os, const TestCase& t)
{
os << t.N;
for (int i : t.phase)
os << " " << i;
os << std::endl;
return os;
}
static const std::vector<TestCase> tests = {
{1, { 0, 1, 2, -1, -1, -1}},
{3, { 0, 1, 2, -1, -1, -1}},
{1, { 0, -1, 1, -1, -1, -1}},
{5, { 0, -1, 1, -1, -1, -1}},
{1, {-1, -1, 0, -1, -1, -1}},
{4, {-1, -1, 0, -1, -1, -1}},
{1, { 0, 1, -1, 2, -1, -1}},
{1, { 0, 1, -1, -1, 2, -1}},
{1, { 0, 1, -1, -1, -1, 2}},
};
#if BOOST_VERSION / 100000 == 1 && BOOST_VERSION / 100 % 1000 > 66
BOOST_DATA_TEST_CASE(RstConvTest, tests)
#else
BOOST_AUTO_TEST_CASE(RstConvTest)
#endif
{
#if BOOST_VERSION / 100000 == 1 && BOOST_VERSION / 100 % 1000 < 67
for (const auto& sample : tests) {
#endif
const auto& cc = Dune::MPIHelper::getCommunication();
Opm::RSTConfig rst;
rst.keywords["CONV"] = sample.N;
std::vector<int> cellMapping(10);
std::iota(cellMapping.begin(), cellMapping.end(), cc.rank()*10);
Dune::BlockVector<Dune::FieldVector<double,6>> residual(10);
// generate data
std::vector<std::vector<int>> max(6);
if (cc.rank() == 0) {
std::random_device rng_device;
std::mt19937 mersenne_engine{rng_device()};
std::uniform_int_distribution<int> dist{0, 10*cc.size()-1};
for (int c = 0; c < 6; ++c) {
if (sample.phase[c] == -1) {
continue;
}
std::vector<int>& v = max[c];
while (v.size() < sample.N) {
int m = dist(mersenne_engine);
if (std::find(v.begin(), v.end(), m) == v.end()) {
v.push_back(m);
}
}
}
}
for (int c = 0; c < 6; ++c) {
std::size_t size = max[c].size();
cc.broadcast(&size, 1, 0);
if (cc.rank() != 0) {
max[c].resize(size);
}
cc.broadcast(max[c].data(), max[c].size(), 0);
}
for (int i = 0; i < 10; ++i) {
for (int c = 0; c < 6; ++c) {
if (sample.phase[c] != -1) {
bool inMax = std::find(max[c].begin(),
max[c].end(),
cellMapping[i]) != max[c].end();
residual[i][sample.phase[c]] = inMax ? 1.0 : 1.0 / (i+2);
}
}
}
Opm::RSTConv cnv([&cellMapping](const int idx) { return cellMapping[idx]; }, cc);
cnv.init(10*cc.size(), rst, sample.phase);
cnv.update(residual);
cnv.update(residual);
if (cc.rank() == 0) {
BOOST_CHECK_EQUAL(cnv.getData().size(), 6);
for (std::size_t i = 0; i < 6; ++i) {
BOOST_CHECK_EQUAL(cnv.getData()[i].size(),
sample.phase[i] == -1 ? 0 : cc.size() * 10);
}
for (int i = 0; i < cc.size() * 10; ++i) {
for (int c = 0; c < 6; ++c) {
if (sample.phase[c] != -1) {
bool inMax = std::find(max[c].begin(),
max[c].end(), i) != max[c].end();
BOOST_CHECK_EQUAL(cnv.getData()[c][i], inMax ? 2 : 0);
}
}
}
}
{
const std::vector<int> conv_new(10, 0);
cnv.updateNewton(conv_new);
const std::vector<int> conv_new_ite1(10, 1);
cnv.updateNewton(conv_new_ite1);
if (cc.rank() == 0) {
BOOST_CHECK_EQUAL(cnv.getConvNew().size(), 10*cc.size());
for (const auto& count : cnv.getConvNew()) {
BOOST_CHECK_EQUAL(count, 2);
}
}
}
#if BOOST_VERSION / 100000 == 1 && BOOST_VERSION / 100 % 1000 < 67
}
#endif
}
int main(int argc, char** argv)
{
Dune::MPIHelper::instance(argc, argv);
#if HAVE_MPI
// register a throwing error handler to allow for
// debugging with "catch throw" in gdb
MPI_Errhandler handler;
MPI_Comm_create_errhandler(MPI_err_handler, &handler);
MPI_Comm_set_errhandler(MPI_COMM_WORLD, handler);
#endif
return boost::unit_test::unit_test_main(&init_unit_test_func, argc, argv);
}
|