BCRSMatrixBlockAssembler.hpp
Go to the documentation of this file.
1/*===========================================================================
2//
3// File: Dune::BCRSMatrixBlockAssembler.hpp
4//
5// Created: 2011-10-02 14:04:38+0200
6//
7// Authors: Ingeborg S. Ligaarden <Ingeborg.Ligaarden@sintef.no>
8// Jostein R. Natvig <Jostein.R.Natvig@sintef.no>
9// Halvor M. Nilsen <HalvorMoll.Nilsen@sintef.no>
10// Atgeirr F. Rasmussen <atgeirr@sintef.no>
11// Bård Skaflestad <Bard.Skaflestad@sintef.no>
12//
13//==========================================================================*/
14
15
16/*
17 Copyright 2011 SINTEF ICT, Applied Mathematics.
18 Copyright 2011 Statoil ASA.
19
20 This file is part of the Open Porous Media Project (OPM).
21
22 OPM is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OPM is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OPM. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPM_BCRSMATRIXBLOCKASSEMBLER_HPP_HEADER
37#define OPM_BCRSMATRIXBLOCKASSEMBLER_HPP_HEADER
38
39#include <cassert>
40#include <cstddef>
41
42#include <dune/common/fvector.hh>
43#include <dune/common/fmatrix.hh>
44#include <dune/istl/bvector.hh>
45#include <dune/istl/bcrsmatrix.hh>
46
47#include <opm/core/transport/implicit/JacobianSystem.hpp>
48
49namespace Opm {
50 namespace ImplicitTransportDefault {
51 namespace ISTLTypeDetails {
52 typedef Dune::FieldVector<double, 1> ScalarVectorBlockType;
53 typedef Dune::FieldMatrix<double, 1, 1> ScalarMatrixBlockType;
54
55 typedef Dune::BlockVector<ScalarVectorBlockType> ScalarBlockVector;
56 typedef Dune::BCRSMatrix <ScalarMatrixBlockType> ScalarBCRSMatrix;
57 }
58
59 template <>
60 class VectorAdder<ISTLTypeDetails::ScalarBlockVector> {
61 public:
62 static void
65 y += x;
66 }
67 };
68
69 template <>
70 class VectorNegater<ISTLTypeDetails::ScalarBlockVector> {
71 public:
72 static void
74 x *= -1.0;
75 }
76 };
77
78 template <>
79 class VectorZero<ISTLTypeDetails::ScalarBlockVector> {
80 public:
81 static void
83 x = 0.0;
84 }
85 };
86
87 template <>
88 class VectorBlockAssembler<ISTLTypeDetails::ScalarBlockVector> {
89 public:
90 template <class Block>
91 static void
92 assemble(::std::size_t ndof,
93 ::std::size_t i ,
94 const Block& b ,
96 assert (ndof == 1); (void) ndof;
97
98 ISTLTypeDetails::ScalarBlockVector::block_type blk(b[0]);
99
100 vec[i] += blk;
101 }
102 };
103
104 template <>
105 class VectorAssign<ISTLTypeDetails::ScalarBlockVector> {
106 public:
107 // y <- x
108 static void
111 y = x;
112 }
113
114 // y <- a*x
115 template <class Scalar>
116 static void
117 assign(const Scalar& a,
120 y = x;
121 y *= a;
122 }
123 };
124
125 template <>
126 class MatrixZero<ISTLTypeDetails::ScalarBCRSMatrix> {
127 public:
128 static void
130 A = 0.0;
131 }
132 };
133
134 template <>
135 class MatrixBlockAssembler<ISTLTypeDetails::ScalarBCRSMatrix>
136 {
137 public:
138 template <class Block>
139 void
140 assembleBlock(size_t ndof, size_t i, size_t j, const Block& b) {
141 assert (ndof == 1); (void) ndof;
142
143 ISTLTypeDetails::ScalarBCRSMatrix::block_type blk(b[0]);
144
145 mat_[i][j] += blk;
146 }
147
148 template <class Connections>
149 void
150 createBlockRow(size_t i, const Connections& conn, size_t ndof) {
151 assert (ndof == 1); (void) ndof;
152 assert (i == i_prev_ + 1); (void) i ;
153
154 ISTLTypeDetails::ScalarBCRSMatrix::CreateIterator ci(mat_, i);
155
156 for (typename Connections::const_iterator
157 c = conn.begin(), e = conn.end(); c != e; ++c) {
158 ci.insert(*c);
159 }
160
161 ++i_prev_;
162 ++ci;
163 }
164
165 void
167
168 void
169 setSize(::std::size_t ndof,
170 ::std::size_t nrow,
171 ::std::size_t ncol,
172 ::std::size_t nnz = 0) {
173
174 assert (ndof == 1); (void) ndof;
175 (void) nnz;
176
177 mat_.setSize (nrow, ncol);
178 mat_.setBuildMode(ISTLTypeDetails::ScalarBCRSMatrix::row_wise);
179
180 i_prev_ = -1;
181 }
182
184 matrix() const { return mat_; }
185
187 matrix() { return mat_; }
188
189 private:
190 ::std::size_t i_prev_;
192 };
193 }
194}
195
196#endif /* OPM_BCRSMATRIXBLOCKASSEMBLER_HPP_HEADER */
ISTLTypeDetails::ScalarBCRSMatrix & matrix()
Definition: BCRSMatrixBlockAssembler.hpp:187
void setSize(::std::size_t ndof, ::std::size_t nrow, ::std::size_t ncol, ::std::size_t nnz=0)
Definition: BCRSMatrixBlockAssembler.hpp:169
void createBlockRow(size_t i, const Connections &conn, size_t ndof)
Definition: BCRSMatrixBlockAssembler.hpp:150
void assembleBlock(size_t ndof, size_t i, size_t j, const Block &b)
Definition: BCRSMatrixBlockAssembler.hpp:140
const ISTLTypeDetails::ScalarBCRSMatrix & matrix() const
Definition: BCRSMatrixBlockAssembler.hpp:184
static void zero(ISTLTypeDetails::ScalarBCRSMatrix &A)
Definition: BCRSMatrixBlockAssembler.hpp:129
static void add(const ISTLTypeDetails::ScalarBlockVector &x, ISTLTypeDetails::ScalarBlockVector &y)
Definition: BCRSMatrixBlockAssembler.hpp:63
static void assign(const Scalar &a, const ISTLTypeDetails::ScalarBlockVector &x, ISTLTypeDetails::ScalarBlockVector &y)
Definition: BCRSMatrixBlockAssembler.hpp:117
static void assign(const ISTLTypeDetails::ScalarBlockVector &x, ISTLTypeDetails::ScalarBlockVector &y)
Definition: BCRSMatrixBlockAssembler.hpp:109
static void assemble(::std::size_t ndof, ::std::size_t i, const Block &b, ISTLTypeDetails::ScalarBlockVector &vec)
Definition: BCRSMatrixBlockAssembler.hpp:92
static void negate(ISTLTypeDetails::ScalarBlockVector &x)
Definition: BCRSMatrixBlockAssembler.hpp:73
static void zero(ISTLTypeDetails::ScalarBlockVector &x)
Definition: BCRSMatrixBlockAssembler.hpp:82
Dune::BlockVector< ScalarVectorBlockType > ScalarBlockVector
Definition: BCRSMatrixBlockAssembler.hpp:55
Dune::FieldMatrix< double, 1, 1 > ScalarMatrixBlockType
Definition: BCRSMatrixBlockAssembler.hpp:53
Dune::BCRSMatrix< ScalarMatrixBlockType > ScalarBCRSMatrix
Definition: BCRSMatrixBlockAssembler.hpp:56
Dune::FieldVector< double, 1 > ScalarVectorBlockType
Definition: BCRSMatrixBlockAssembler.hpp:52
Definition: BlackoilFluid.hpp:32