3 #ifndef DUNE_ISTL_BTDMATRIX_HH
4 #define DUNE_ISTL_BTDMATRIX_HH
6 #include <dune/common/fmatrix.hh>
7 #include <dune/common/scalarvectorview.hh>
8 #include <dune/common/scalarmatrixview.hh>
27 template <
class B,
class A=std::allocator<B> >
35 using field_type =
typename Imp::BlockTraits<B>::field_type;
50 [[deprecated(
"Use free blockLevel function. Will be removed after 2.8.")]]
61 for (
size_t i=0; i<size; i++)
67 for (
size_t i=0; i<size; i++) {
81 auto nonZeros = (size==0) ? 0 : size + 2*(size-1);
88 for (
size_t i=0; i<size; i++)
94 for (
size_t i=0; i<size; i++) {
123 void solve (V& x,
const V& rhs)
const {
127 auto&& x0 = Impl::asVector(x[0]);
128 auto&& rhs0 = Impl::asVector(rhs[0]);
129 Impl::asMatrix((*
this)[0][0]).solve(x0, rhs0);
135 std::vector<block_type> c(this->
N()-1);
136 for (
size_t i=0; i<this->
N()-1; i++)
137 c[i] = (*
this)[i][i+1];
141 Impl::asMatrix(a_00_inv).invert();
145 Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[0]));
150 auto&& d_0 = Impl::asVector(d[0]);
151 Impl::asMatrix(a_00_inv).mv(Impl::asVector(d_0_tmp),d_0);
153 for (
unsigned int i = 1; i < this->
N(); i++) {
157 tmp = (*this)[i][i-1];
158 Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[i-1]));
162 Impl::asMatrix(
id).invert();
165 Impl::asMatrix(c[i]).leftmultiply(Impl::asMatrix(
id));
169 auto&& d_i = Impl::asVector(d[i]);
170 Impl::asMatrix((*
this)[i][i-1]).mmv(Impl::asVector(d[i-1]), d_i);
172 Impl::asMatrix(
id).mv(Impl::asVector(tmpVec), d_i);
176 x[this->
N() - 1] = d[this->
N() - 1];
177 for (
int i = this->
N() - 2; i >= 0; i--) {
180 auto&& x_i = Impl::asVector(x[i]);
181 Impl::asMatrix(c[i]).mmv(Impl::asVector(x[i+1]), x_i);
197 void endrowsizes () {}
199 void endindices () {}
Implementation of the BCRSMatrix class.
Helper functions for determining the vector/matrix block level.
Col col
Definition: matrixmatrix.hh:349
Definition: allocator.hh:9
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:464
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:909
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1147
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:528
void addindex(size_type row, size_type col)
add index (row,col) to the matrix
Definition: bcrsmatrix.hh:1189
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1246
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1970
void setSize(size_type rows, size_type columns, size_type nnz=0)
Set the size of the matrix.
Definition: bcrsmatrix.hh:859
A block-tridiagonal matrix.
Definition: btdmatrix.hh:29
static constexpr auto blocklevel
increment block level counter
Definition: btdmatrix.hh:51
BTDMatrix(size_type size)
Definition: btdmatrix.hh:56
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:123
BTDMatrix & operator=(const BTDMatrix &other)
assignment
Definition: btdmatrix.hh:106
A::size_type size_type
implement row_type with compressed vector
Definition: btdmatrix.hh:47
A allocator_type
export the allocator type
Definition: btdmatrix.hh:41
B block_type
export the type representing the components
Definition: btdmatrix.hh:38
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: btdmatrix.hh:35
BTDMatrix()
Default constructor.
Definition: btdmatrix.hh:54
void setSize(size_type size)
Resize the matrix. Invalidates the content!
Definition: btdmatrix.hh:79