JacobianSystem.hpp
Go to the documentation of this file.
1/*===========================================================================
2//
3// File: JacobianSystem.hpp
4//
5// Created: 2011-09-30 19:23:31+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_JACOBIANSYSTEM_HPP_HEADER
37#define OPM_JACOBIANSYSTEM_HPP_HEADER
38
39#include <cassert>
40#include <cmath>
41#include <cstddef>
42
43#include <algorithm>
44//#include <array>
45#include <functional>
46#include <numeric>
47
48namespace Opm {
49 namespace ImplicitTransportDefault {
50 template <class BaseVec>
52 public:
53 // y += x
54 static void
55 add(const BaseVec& x, BaseVec& y) {
56 typedef typename BaseVec::value_type VT;
57 ::std::ranges::transform(x, y, y.begin(), ::std::plus<VT>());
58 }
59 };
60
61 template <class BaseVec>
63 public:
64 // x *= -1
65 static void
66 negate(BaseVec& x) {
67 typedef typename BaseVec::value_type VT;
68 ::std::ranges::transform(x, x.begin(), ::std::negate<VT>());
69 }
70 };
71
72 template <class BaseVec>
73 class VectorZero {
74 public:
75 static void
76 zero(BaseVec& x) {
77 typedef typename BaseVec::value_type VT;
78
79 ::std::ranges::fill(x, VT(0.0));
80 }
81 };
82
83 template <class BaseVec>
85 public:
86 // y <- x
87 static void
88 assign(const BaseVec& x, BaseVec& y) {
89 ::std::ranges::copy(x, y.begin());
90 }
91
92 // y <- a*x
93 template <class Scalar>
94 static void
95 assign(const Scalar& a, const BaseVec& x, BaseVec& y) {
96 ::std::ranges::transform(x, y.begin(),
97 [&a](auto z) { return a * z; });
98 }
99 };
100
101 template <class Matrix>
103
104 template <class BaseVec>
106 public:
107 template <class Block>
108 static void
109 assemble(::std::size_t ndof,
110 ::std::size_t i ,
111 const Block& b ,
112 BaseVec& vec ) {
113
114 for (::std::size_t d = 0; d < ndof; ++d) {
115 vec[i*ndof + d] += b[d];
116 }
117 }
118 };
119
120 template <class BaseVec>
122 public:
123 explicit VectorSizeSetter(BaseVec& v) : v_(v) {}
124
125 void
126 setSize(::std::size_t ndof, ::std::size_t m) {
127 v_.resize(ndof * m);
128 }
129
130 private:
131 BaseVec& v_;
132 };
133
134 template <class BaseVec ,
135 template <class> class VSzSetter = VectorSizeSetter ,
136 template <class> class VAdd = VectorAdder ,
137 template <class> class VBlkAsm = VectorBlockAssembler>
139 enum { Residual = 0, Increment = 1, Solution = 2 };
140
141 public:
142 void
143 setSize(::std::size_t ndof, ::std::size_t m) {
144#if 0
145 typedef typename ::std::array<BaseVec, 3>::iterator VAI;
146
147 for (VAI i = vcoll_.begin(), e = vcoll_.end(); i != e; ++i) {
148 VSzSetter<BaseVec>(*i).setSize(ndof, m);
149 }
150#endif
151 for (::std::size_t i = 0; i < sizeof (vcoll_) / sizeof (vcoll_[0]); ++i) {
152 VSzSetter<BaseVec>(vcoll_[i]).setSize(ndof, m);
153 }
154
155 ndof_ = ndof;
156 }
157
158 void
160 VAdd<BaseVec>::add(vcoll_[ Increment ], vcoll_[ Solution ]);
161 }
162
163 template <class Block>
164 void
165 assembleBlock(::std::size_t ndof,
166 ::std::size_t i ,
167 const Block& b ) {
168
169 assert (ndof_ > 0 );
170 assert (ndof == ndof_);
171
172 VBlkAsm<BaseVec>::assemble(ndof, i, b, vcoll_[ Residual ]);
173 }
174
175 typedef BaseVec vector_type;
176
177 const vector_type& increment() const { return vcoll_[ Increment ]; }
178 const vector_type& residual () const { return vcoll_[ Residual ]; }
179 const vector_type& solution () const { return vcoll_[ Solution ]; }
180
181 // Write access for Newton solver purposes
182 vector_type& writableIncrement() { return vcoll_[ Increment ]; }
183 vector_type& writableResidual () { return vcoll_[ Residual ]; }
184 vector_type& writableSolution () { return vcoll_[ Solution ]; }
185
186 private:
187 ::std::size_t ndof_ ;
188 BaseVec vcoll_[3];
189 //::std::array<BaseVec, 3> vcoll_;
190 };
191
192 template <class Matrix>
193 class MatrixBlockAssembler
194 /* {
195 public:
196 template <class Block>
197 void
198 assembleBlock(size_t n, size_t i, size j, const Block& b);
199
200 template <class Connections>
201 void
202 createBlockRow(size_t i, const Connections& conn, size_t n);
203
204 void
205 finalizeStructure();
206
207 void
208 setSize(size_t ndof, size_t m, size_t n, size_t nnz = 0);
209
210 Matrix& matrix() ;
211 const Matrix& matrix() const;
212 } */;
213
214 template <class Matrix ,
215 class NVecCollection>
217 public:
219
221 typedef MatrixBlockAssembler<Matrix> assembler_type;
222 typedef NVecCollection vector_collection_type;
223 typedef typename NVecCollection::vector_type vector_type;
224
225 assembler_type& matasm() { return mba_ ; }
226 NVecCollection& vector() { return sysvec_ ; }
227 const matrix_type& matrix() const { return mba_.matrix(); }
228 matrix_type& writableMatrix() { return mba_.matrix(); }
229
230 void
231 setSize(::std::size_t ndof,
232 ::std::size_t m ,
233 ::std::size_t nnz = 0) {
234
235 mba_ .setSize(ndof, m, m, nnz);
236 sysvec_.setSize(ndof, m );
237 }
238
239 private:
241 JacobianSystem& operator=(const JacobianSystem&);
242
243 MatrixBlockAssembler<Matrix> mba_ ; // Coefficient matrix
244 NVecCollection sysvec_; // Residual
245 };
246 }
247}
248
249#endif /* OPM_JACOBIANSYSTEM_HPP_HEADER */
Definition: JacobianSystem.hpp:216
Matrix matrix_type
Definition: JacobianSystem.hpp:220
JacobianSystem()
Definition: JacobianSystem.hpp:218
NVecCollection vector_collection_type
Definition: JacobianSystem.hpp:222
matrix_type & writableMatrix()
Definition: JacobianSystem.hpp:228
void setSize(::std::size_t ndof, ::std::size_t m, ::std::size_t nnz=0)
Definition: JacobianSystem.hpp:231
NVecCollection & vector()
Definition: JacobianSystem.hpp:226
MatrixBlockAssembler< Matrix > assembler_type
Definition: JacobianSystem.hpp:221
NVecCollection::vector_type vector_type
Definition: JacobianSystem.hpp:223
assembler_type & matasm()
Definition: JacobianSystem.hpp:225
const matrix_type & matrix() const
Definition: JacobianSystem.hpp:227
Definition: JacobianSystem.hpp:102
void assembleBlock(::std::size_t ndof, ::std::size_t i, const Block &b)
Definition: JacobianSystem.hpp:165
void setSize(::std::size_t ndof, ::std::size_t m)
Definition: JacobianSystem.hpp:143
const vector_type & increment() const
Definition: JacobianSystem.hpp:177
vector_type & writableSolution()
Definition: JacobianSystem.hpp:184
void addIncrement()
Definition: JacobianSystem.hpp:159
vector_type & writableIncrement()
Definition: JacobianSystem.hpp:182
vector_type & writableResidual()
Definition: JacobianSystem.hpp:183
const vector_type & solution() const
Definition: JacobianSystem.hpp:179
BaseVec vector_type
Definition: JacobianSystem.hpp:175
const vector_type & residual() const
Definition: JacobianSystem.hpp:178
Definition: JacobianSystem.hpp:51
static void add(const BaseVec &x, BaseVec &y)
Definition: JacobianSystem.hpp:55
Definition: JacobianSystem.hpp:84
static void assign(const Scalar &a, const BaseVec &x, BaseVec &y)
Definition: JacobianSystem.hpp:95
static void assign(const BaseVec &x, BaseVec &y)
Definition: JacobianSystem.hpp:88
Definition: JacobianSystem.hpp:105
static void assemble(::std::size_t ndof, ::std::size_t i, const Block &b, BaseVec &vec)
Definition: JacobianSystem.hpp:109
Definition: JacobianSystem.hpp:62
static void negate(BaseVec &x)
Definition: JacobianSystem.hpp:66
Definition: JacobianSystem.hpp:121
VectorSizeSetter(BaseVec &v)
Definition: JacobianSystem.hpp:123
void setSize(::std::size_t ndof, ::std::size_t m)
Definition: JacobianSystem.hpp:126
Definition: JacobianSystem.hpp:73
static void zero(BaseVec &x)
Definition: JacobianSystem.hpp:76
Dune::BCRSMatrix< Dune::FieldMatrix< double, 1, 1 > > Matrix
A sparse matrix holding our operator.
Definition: matrixops.hpp:27
Definition: ImplicitAssembly.hpp:43