Dune Core Modules (2.9.0)

btdmatrix.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_ISTL_BTDMATRIX_HH
6 #define DUNE_ISTL_BTDMATRIX_HH
7 
8 #include <dune/common/fmatrix.hh>
11 #include <dune/istl/bcrsmatrix.hh>
12 #include <dune/istl/blocklevel.hh>
13 
19 namespace Dune {
29  template <class B, class A=std::allocator<B> >
30  class BTDMatrix : public BCRSMatrix<B,A>
31  {
32  public:
33 
34  //===== type definitions and constants
35 
37  using field_type = typename Imp::BlockTraits<B>::field_type;
38 
40  typedef B block_type;
41 
43  typedef A allocator_type;
44 
46  //typedef BCRSMatrix<B,A>::row_type row_type;
47 
49  typedef typename A::size_type size_type;
50 
52  [[deprecated("Use free blockLevel function. Will be removed after 2.8.")]]
53  static constexpr auto blocklevel = blockLevel<B>()+1;
54 
56  BTDMatrix() : BCRSMatrix<B,A>() {}
57 
58  explicit BTDMatrix(size_type size)
59  : BCRSMatrix<B,A>(size, size, BCRSMatrix<B,A>::random)
60  {
61  // Set number of entries for each row
62  // All rows get three entries, except for the first and the last one
63  for (size_t i=0; i<size; i++)
64  this->BCRSMatrix<B,A>::setrowsize(i, 3 - (i==0) - (i==(size-1)));
65 
67 
68  // The actual entries for each row
69  for (size_t i=0; i<size; i++) {
70  if (i>0)
71  this->BCRSMatrix<B,A>::addindex(i, i-1);
72  this->BCRSMatrix<B,A>::addindex(i, i );
73  if (i<size-1)
74  this->BCRSMatrix<B,A>::addindex(i, i+1);
75  }
76 
78  }
79 
81  void setSize(size_type size)
82  {
83  auto nonZeros = (size==0) ? 0 : size + 2*(size-1);
84  this->BCRSMatrix<B,A>::setSize(size, // rows
85  size, // columns
86  nonZeros);
87 
88  // Set number of entries for each row
89  // All rows get three entries, except for the first and the last one
90  for (size_t i=0; i<size; i++)
91  this->BCRSMatrix<B,A>::setrowsize(i, 3 - (i==0) - (i==(size-1)));
92 
94 
95  // The actual entries for each row
96  for (size_t i=0; i<size; i++) {
97  if (i>0)
98  this->BCRSMatrix<B,A>::addindex(i, i-1);
99  this->BCRSMatrix<B,A>::addindex(i, i );
100  if (i<size-1)
101  this->BCRSMatrix<B,A>::addindex(i, i+1);
102  }
103 
105  }
106 
108  BTDMatrix& operator= (const BTDMatrix& other) {
109  this->BCRSMatrix<B,A>::operator=(other);
110  return *this;
111  }
112 
116  return *this;
117  }
118 
124  template <class V>
125  void solve (V& x, const V& rhs) const {
126 
127  // special handling for 1x1 matrices. The generic algorithm doesn't work for them
128  if (this->N()==1) {
129  auto&& x0 = Impl::asVector(x[0]);
130  auto&& rhs0 = Impl::asVector(rhs[0]);
131  Impl::asMatrix((*this)[0][0]).solve(x0, rhs0);
132  return;
133  }
134 
135  // Make copies of the rhs and the right matrix band
136  V d = rhs;
137  std::vector<block_type> c(this->N()-1);
138  for (size_t i=0; i<this->N()-1; i++)
139  c[i] = (*this)[i][i+1];
140 
141  /* Modify the coefficients. */
142  block_type a_00_inv = (*this)[0][0];
143  Impl::asMatrix(a_00_inv).invert();
144 
145  //c[0] /= (*this)[0][0]; /* Division by zero risk. */
146  block_type tmp = a_00_inv;
147  Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[0]));
148  c[0] = tmp;
149 
150  // d = a^{-1} d /* Division by zero would imply a singular matrix. */
151  auto d_0_tmp = d[0];
152  auto&& d_0 = Impl::asVector(d[0]);
153  Impl::asMatrix(a_00_inv).mv(Impl::asVector(d_0_tmp),d_0);
154 
155  for (unsigned int i = 1; i < this->N(); i++) {
156 
157  // id = ( a_ii - c_{i-1} a_{i, i-1} ) ^{-1}
158  block_type tmp;
159  tmp = (*this)[i][i-1];
160  Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[i-1]));
161 
162  block_type id = (*this)[i][i];
163  id -= tmp;
164  Impl::asMatrix(id).invert(); /* Division by zero risk. */
165 
166  if (i<c.size()) {
167  Impl::asMatrix(c[i]).leftmultiply(Impl::asMatrix(id)); /* Last value calculated is redundant. */
168  }
169 
170  // d[i] = (d[i] - d[i-1] * (*this)[i][i-1]) * id;
171  auto&& d_i = Impl::asVector(d[i]);
172  Impl::asMatrix((*this)[i][i-1]).mmv(Impl::asVector(d[i-1]), d_i);
173  auto tmpVec = d[i];
174  Impl::asMatrix(id).mv(Impl::asVector(tmpVec), d_i);
175  }
176 
177  /* Now back substitute. */
178  x[this->N() - 1] = d[this->N() - 1];
179  for (int i = this->N() - 2; i >= 0; i--) {
180  //x[i] = d[i] - c[i] * x[i + 1];
181  x[i] = d[i];
182  auto&& x_i = Impl::asVector(x[i]);
183  Impl::asMatrix(c[i]).mmv(Impl::asVector(x[i+1]), x_i);
184  }
185 
186  }
187 
188  private:
189 
190  // ////////////////////////////////////////////////////////////////////////////
191  // The following methods from the base class should now actually be called
192  // ////////////////////////////////////////////////////////////////////////////
193 
194  // createbegin and createend should be in there, too, but I can't get it to compile
195  // BCRSMatrix<B,A>::CreateIterator createbegin () {}
196  // BCRSMatrix<B,A>::CreateIterator createend () {}
197  void setrowsize (size_type i, size_type s) {}
198  void incrementrowsize (size_type i) {}
199  void endrowsizes () {}
200  void addindex (size_type row, size_type col) {}
201  void endindices () {}
202  };
203 
204  template<typename B, typename A>
205  struct FieldTraits< BTDMatrix<B, A> >
206  {
207  using field_type = typename BTDMatrix<B, A>::field_type;
208  using real_type = typename FieldTraits<field_type>::real_type;
209  };
210 
213 } // end namespace Dune
214 
215 #endif
Implementation of the BCRSMatrix class.
Helper functions for determining the vector/matrix block level.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:911
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1149
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:530
void addindex(size_type row, size_type col)
add index (row,col) to the matrix
Definition: bcrsmatrix.hh:1191
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1248
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1972
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:861
A block-tridiagonal matrix.
Definition: btdmatrix.hh:31
static constexpr auto blocklevel
increment block level counter
Definition: btdmatrix.hh:53
void solve(V &x, const V &rhs) const
Use the Thomas algorithm to solve the system Ax=b in O(n) time.
Definition: btdmatrix.hh:125
BTDMatrix & operator=(const BTDMatrix &other)
assignment
Definition: btdmatrix.hh:108
A::size_type size_type
implement row_type with compressed vector
Definition: btdmatrix.hh:49
A allocator_type
export the allocator type
Definition: btdmatrix.hh:43
B block_type
export the type representing the components
Definition: btdmatrix.hh:40
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: btdmatrix.hh:37
BTDMatrix()
Default constructor.
Definition: btdmatrix.hh:56
void setSize(size_type size)
Resize the matrix. Invalidates the content!
Definition: btdmatrix.hh:81
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Dune namespace.
Definition: alignedallocator.hh:13
Implements a scalar matrix view wrapper around an existing scalar.
Implements a scalar vector view wrapper around an existing scalar.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 27, 22:29, 2024)