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
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
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 2 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/>.
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
#ifndef EWOMS_MATRIX_BLOCK_HH
#define EWOMS_MATRIX_BLOCK_HH
#include <dune/common/dynmatrix.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/version.hh>
#include <dune/common/typetraits.hh>
#include <dune/istl/superlu.hh>
#include <dune/istl/umfpack.hh>
#include <dune/istl/matrixutils.hh>
#include <opm/common/Exceptions.hpp>
#include <limits>
namespace Opm {
namespace detail {
template <typename K, int m, int n>
static inline void invertMatrix(Dune::FieldMatrix<K,m,n>& matrix)
{
matrix.invert();
}
template <typename K>
static inline void invertMatrix(Dune::FieldMatrix<K,1,1>& matrix)
{
Dune::FieldMatrix<K,1,1> tmp(matrix);
Dune::FMatrixHelp::invertMatrix(tmp,matrix);
}
template <typename K>
static inline void invertMatrix(Dune::FieldMatrix<K,2,2>& matrix)
{
Dune::FieldMatrix<K,2,2> tmp(matrix);
Dune::FMatrixHelp::invertMatrix(tmp,matrix);
}
template <typename K>
static inline void invertMatrix(Dune::FieldMatrix<K,3,3>& matrix)
{
Dune::FieldMatrix<K,3,3> tmp(matrix);
Dune::FMatrixHelp::invertMatrix(tmp,matrix);
}
//! invert 4x4 Matrix without changing the original matrix
template <template<class K> class Matrix, typename K>
static inline K invertMatrix4(const Matrix<K>& matrix, Matrix<K>& inverse)
{
inverse[0][0] = matrix[1][1] * matrix[2][2] * matrix[3][3] -
matrix[1][1] * matrix[2][3] * matrix[3][2] -
matrix[2][1] * matrix[1][2] * matrix[3][3] +
matrix[2][1] * matrix[1][3] * matrix[3][2] +
matrix[3][1] * matrix[1][2] * matrix[2][3] -
matrix[3][1] * matrix[1][3] * matrix[2][2];
inverse[1][0] = -matrix[1][0] * matrix[2][2] * matrix[3][3] +
matrix[1][0] * matrix[2][3] * matrix[3][2] +
matrix[2][0] * matrix[1][2] * matrix[3][3] -
matrix[2][0] * matrix[1][3] * matrix[3][2] -
matrix[3][0] * matrix[1][2] * matrix[2][3] +
matrix[3][0] * matrix[1][3] * matrix[2][2];
inverse[2][0] = matrix[1][0] * matrix[2][1] * matrix[3][3] -
matrix[1][0] * matrix[2][3] * matrix[3][1] -
matrix[2][0] * matrix[1][1] * matrix[3][3] +
matrix[2][0] * matrix[1][3] * matrix[3][1] +
matrix[3][0] * matrix[1][1] * matrix[2][3] -
matrix[3][0] * matrix[1][3] * matrix[2][1];
inverse[3][0] = -matrix[1][0] * matrix[2][1] * matrix[3][2] +
matrix[1][0] * matrix[2][2] * matrix[3][1] +
matrix[2][0] * matrix[1][1] * matrix[3][2] -
matrix[2][0] * matrix[1][2] * matrix[3][1] -
matrix[3][0] * matrix[1][1] * matrix[2][2] +
matrix[3][0] * matrix[1][2] * matrix[2][1];
inverse[0][1]= -matrix[0][1] * matrix[2][2] * matrix[3][3] +
matrix[0][1] * matrix[2][3] * matrix[3][2] +
matrix[2][1] * matrix[0][2] * matrix[3][3] -
matrix[2][1] * matrix[0][3] * matrix[3][2] -
matrix[3][1] * matrix[0][2] * matrix[2][3] +
matrix[3][1] * matrix[0][3] * matrix[2][2];
inverse[1][1] = matrix[0][0] * matrix[2][2] * matrix[3][3] -
matrix[0][0] * matrix[2][3] * matrix[3][2] -
matrix[2][0] * matrix[0][2] * matrix[3][3] +
matrix[2][0] * matrix[0][3] * matrix[3][2] +
matrix[3][0] * matrix[0][2] * matrix[2][3] -
matrix[3][0] * matrix[0][3] * matrix[2][2];
inverse[2][1] = -matrix[0][0] * matrix[2][1] * matrix[3][3] +
matrix[0][0] * matrix[2][3] * matrix[3][1] +
matrix[2][0] * matrix[0][1] * matrix[3][3] -
matrix[2][0] * matrix[0][3] * matrix[3][1] -
matrix[3][0] * matrix[0][1] * matrix[2][3] +
matrix[3][0] * matrix[0][3] * matrix[2][1];
inverse[3][1] = matrix[0][0] * matrix[2][1] * matrix[3][2] -
matrix[0][0] * matrix[2][2] * matrix[3][1] -
matrix[2][0] * matrix[0][1] * matrix[3][2] +
matrix[2][0] * matrix[0][2] * matrix[3][1] +
matrix[3][0] * matrix[0][1] * matrix[2][2] -
matrix[3][0] * matrix[0][2] * matrix[2][1];
inverse[0][2] = matrix[0][1] * matrix[1][2] * matrix[3][3] -
matrix[0][1] * matrix[1][3] * matrix[3][2] -
matrix[1][1] * matrix[0][2] * matrix[3][3] +
matrix[1][1] * matrix[0][3] * matrix[3][2] +
matrix[3][1] * matrix[0][2] * matrix[1][3] -
matrix[3][1] * matrix[0][3] * matrix[1][2];
inverse[1][2] = -matrix[0][0] * matrix[1][2] * matrix[3][3] +
matrix[0][0] * matrix[1][3] * matrix[3][2] +
matrix[1][0] * matrix[0][2] * matrix[3][3] -
matrix[1][0] * matrix[0][3] * matrix[3][2] -
matrix[3][0] * matrix[0][2] * matrix[1][3] +
matrix[3][0] * matrix[0][3] * matrix[1][2];
inverse[2][2] = matrix[0][0] * matrix[1][1] * matrix[3][3] -
matrix[0][0] * matrix[1][3] * matrix[3][1] -
matrix[1][0] * matrix[0][1] * matrix[3][3] +
matrix[1][0] * matrix[0][3] * matrix[3][1] +
matrix[3][0] * matrix[0][1] * matrix[1][3] -
matrix[3][0] * matrix[0][3] * matrix[1][1];
inverse[3][2] = -matrix[0][0] * matrix[1][1] * matrix[3][2] +
matrix[0][0] * matrix[1][2] * matrix[3][1] +
matrix[1][0] * matrix[0][1] * matrix[3][2] -
matrix[1][0] * matrix[0][2] * matrix[3][1] -
matrix[3][0] * matrix[0][1] * matrix[1][2] +
matrix[3][0] * matrix[0][2] * matrix[1][1];
inverse[0][3] = -matrix[0][1] * matrix[1][2] * matrix[2][3] +
matrix[0][1] * matrix[1][3] * matrix[2][2] +
matrix[1][1] * matrix[0][2] * matrix[2][3] -
matrix[1][1] * matrix[0][3] * matrix[2][2] -
matrix[2][1] * matrix[0][2] * matrix[1][3] +
matrix[2][1] * matrix[0][3] * matrix[1][2];
inverse[1][3] = matrix[0][0] * matrix[1][2] * matrix[2][3] -
matrix[0][0] * matrix[1][3] * matrix[2][2] -
matrix[1][0] * matrix[0][2] * matrix[2][3] +
matrix[1][0] * matrix[0][3] * matrix[2][2] +
matrix[2][0] * matrix[0][2] * matrix[1][3] -
matrix[2][0] * matrix[0][3] * matrix[1][2];
inverse[2][3] = -matrix[0][0] * matrix[1][1] * matrix[2][3] +
matrix[0][0] * matrix[1][3] * matrix[2][1] +
matrix[1][0] * matrix[0][1] * matrix[2][3] -
matrix[1][0] * matrix[0][3] * matrix[2][1] -
matrix[2][0] * matrix[0][1] * matrix[1][3] +
matrix[2][0] * matrix[0][3] * matrix[1][1];
inverse[3][3] = matrix[0][0] * matrix[1][1] * matrix[2][2] -
matrix[0][0] * matrix[1][2] * matrix[2][1] -
matrix[1][0] * matrix[0][1] * matrix[2][2] +
matrix[1][0] * matrix[0][2] * matrix[2][1] +
matrix[2][0] * matrix[0][1] * matrix[1][2] -
matrix[2][0] * matrix[0][2] * matrix[1][1];
K det = matrix[0][0] * inverse[0][0] + matrix[0][1] * inverse[1][0] +
matrix[0][2] * inverse[2][0] + matrix[0][3] * inverse[3][0];
// return identity for singular or nearly singular matrices.
if (std::abs(det) < 1e-40) {
inverse = std::numeric_limits<K>::quiet_NaN();
throw NumericalProblem("Singular matrix");
} else
inverse *= 1.0 / det;
return det;
}
template<class K> using FMat4 = Dune::FieldMatrix<K,4,4>;
template <typename K>
static inline void invertMatrix(Dune::FieldMatrix<K,4,4>& matrix)
{
FMat4<K> tmp(matrix);
invertMatrix4<FMat4>(tmp, matrix);
}
template <typename K>
static inline void invertMatrix(Dune::DynamicMatrix<K>& matrix)
{
// this function is only for 4 X 4 matrix
// for 4 X 4 matrix, using the invertMatrix() function above
// it is for temporary usage, mainly to reduce the huge burden of testing
// what algorithm should be used to invert 4 X 4 matrix will be handled
// as a seperate issue
if (matrix.rows() == 4) {
Dune::DynamicMatrix<K> A = matrix;
invertMatrix4(A, matrix);
return;
}
#if ! DUNE_VERSION_NEWER( DUNE_COMMON, 2, 7 )
Dune::FMatrixPrecision<K>::set_singular_limit(1.e-30);
#endif
matrix.invert();
}
} // namespace detail
template <class Scalar, int n, int m>
class MatrixBlock : public Dune::FieldMatrix<Scalar, n, m>
{
public:
using BaseType = Dune::FieldMatrix<Scalar, n, m> ;
using BaseType::operator= ;
using BaseType::rows;
using BaseType::cols;
MatrixBlock()
: BaseType(Scalar(0.0))
{}
explicit MatrixBlock(const Scalar value)
: BaseType(value)
{}
void invert()
{ detail::invertMatrix(asBase()); }
const BaseType& asBase() const
{ return static_cast<const BaseType&>(*this); }
BaseType& asBase()
{ return static_cast<BaseType&>(*this); }
};
} // namespace Opm
namespace Dune {
template<class K, int n, int m>
void print_row(std::ostream& s, const Opm::MatrixBlock<K, n, m>& A,
typename FieldMatrix<K, n, m>::size_type I,
typename FieldMatrix<K, n, m>::size_type J,
typename FieldMatrix<K, n, m>::size_type therow,
int width,
int precision)
{ print_row(s, A.asBase(), I, J, therow, width, precision); }
template <typename Scalar, int n, int m>
struct MatrixDimension<Opm::MatrixBlock<Scalar, n, m> >
: public MatrixDimension<typename Opm::MatrixBlock<Scalar, n, m>::BaseType>
{ };
#if HAVE_UMFPACK
/// \brief UMFPack specialization for Opm::MatrixBlock to make AMG happy
///
/// Without this the empty default implementation would be used.
template <typename T, typename A, int n, int m>
class UMFPack<BCRSMatrix<Opm::MatrixBlock<T, n, m>, A> >
: public UMFPack<BCRSMatrix<FieldMatrix<T, n, m>, A> >
{
using Base = UMFPack<BCRSMatrix<FieldMatrix<T, n, m>, A> >;
using Matrix = BCRSMatrix<FieldMatrix<T, n, m>, A>;
public:
using RealMatrix = BCRSMatrix<Opm::MatrixBlock<T, n, m>, A>;
UMFPack(const RealMatrix& matrix, int verbose, bool)
: Base(reinterpret_cast<const Matrix&>(matrix), verbose)
{}
};
#endif
#if HAVE_SUPERLU
/// \brief SuperLU specialization for Opm::MatrixBlock to make AMG happy
///
/// Without this the empty default implementation would be used.
template <typename T, typename A, int n, int m>
class SuperLU<BCRSMatrix<Opm::MatrixBlock<T, n, m>, A> >
: public SuperLU<BCRSMatrix<FieldMatrix<T, n, m>, A> >
{
using Base = SuperLU<BCRSMatrix<FieldMatrix<T, n, m>, A> >;
using Matrix = BCRSMatrix<FieldMatrix<T, n, m>, A>;
public:
using RealMatrix = BCRSMatrix<Opm::MatrixBlock<T, n, m>, A>;
SuperLU(const RealMatrix& matrix, int verb, bool reuse=true)
: Base(reinterpret_cast<const Matrix&>(matrix), verb, reuse)
{}
};
#endif
template<typename T, int n, int m>
struct IsNumber<Opm::MatrixBlock<T, n, m>>
: public IsNumber<Dune::FieldMatrix<T,n,m>>
{};
} // end namespace Dune
#endif
|