5#ifndef DUNE_ISTL_BTDMATRIX_HH
6#define DUNE_ISTL_BTDMATRIX_HH
8#include <dune/common/fmatrix.hh>
9#include <dune/common/scalarvectorview.hh>
10#include <dune/common/scalarmatrixview.hh>
29 template <
class B,
class A=std::allocator<B> >
37 using field_type =
typename Imp::BlockTraits<B>::field_type;
59 for (
size_t i=0; i<size; i++)
65 for (
size_t i=0; i<size; i++) {
79 auto nonZeros = (size==0) ? 0 : size + 2*(size-1);
86 for (
size_t i=0; i<size; i++)
92 for (
size_t i=0; i<size; i++) {
121 void solve (V& x,
const V& rhs)
const {
125 auto&& x0 = Impl::asVector(x[0]);
126 auto&& rhs0 = Impl::asVector(rhs[0]);
127 Impl::asMatrix((*
this)[0][0]).solve(x0, rhs0);
133 std::vector<block_type> c(this->N()-1);
134 for (
size_t i=0; i<this->N()-1; i++)
135 c[i] = (*
this)[i][i+1];
139 Impl::asMatrix(a_00_inv).invert();
143 Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[0]));
148 auto&& d_0 = Impl::asVector(d[0]);
149 Impl::asMatrix(a_00_inv).mv(Impl::asVector(d_0_tmp),d_0);
151 for (
unsigned int i = 1; i < this->N(); i++) {
155 tmp = (*this)[i][i-1];
156 Impl::asMatrix(tmp).rightmultiply(Impl::asMatrix(c[i-1]));
160 Impl::asMatrix(
id).invert();
163 Impl::asMatrix(c[i]).leftmultiply(Impl::asMatrix(
id));
167 auto&& d_i = Impl::asVector(d[i]);
168 Impl::asMatrix((*
this)[i][i-1]).mmv(Impl::asVector(d[i-1]), d_i);
170 Impl::asMatrix(
id).mv(Impl::asVector(tmpVec), d_i);
174 x[this->N() - 1] = d[this->N() - 1];
175 for (
int i = this->N() - 2; i >= 0; i--) {
178 auto&& x_i = Impl::asVector(x[i]);
179 Impl::asMatrix(c[i]).mmv(Impl::asVector(x[i+1]), x_i);
195 void endrowsizes () {}
197 void endindices () {}
200 template<
typename B,
typename A>
204 using real_type =
typename FieldTraits<field_type>::real_type;
Implementation of the BCRSMatrix class.
Helper functions for determining the vector/matrix block level.
*The type for the index access and the size typedef A::size_type size_type
Definition bcrsmatrix.hh:497
@ random
Build entries randomly.
Definition bcrsmatrix.hh:526
Definition allocator.hh:11
MatrixSparsityPatternGatherScatter< M, I >::ColIter MatrixSparsityPatternGatherScatter< M, I >::col
Definition matrixredistribute.hh:591
A block-tridiagonal matrix.
Definition btdmatrix.hh:31
BTDMatrix(size_type size)
Definition btdmatrix.hh:54
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:121
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 & operator=(const BTDMatrix &other)
assignment
Definition btdmatrix.hh:104
BTDMatrix()
Default constructor.
Definition btdmatrix.hh:52
void setSize(size_type size)
Resize the matrix. Invalidates the content!
Definition btdmatrix.hh:77
typename FieldTraits< field_type >::real_type real_type
Definition btdmatrix.hh:204
typename BTDMatrix< B, A >::field_type field_type
Definition btdmatrix.hh:203
Definition matrixutils.hh:24