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