DUNE PDELab (2.8)

bdmatrix.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ISTL_BDMATRIX_HH
4#define DUNE_ISTL_BDMATRIX_HH
5
6#include <memory>
7
10
13
19namespace Dune {
29 template <class B, class A=std::allocator<B> >
30 class BDMatrix : 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 function blockLevel(). Will be removed after 2.8.")]]
53 static constexpr unsigned int blocklevel = blockLevel<B>()+1;
54
56 BDMatrix() : BCRSMatrix<B,A>() {}
57
58 explicit BDMatrix(int size)
59 : BCRSMatrix<B,A>(size, size, BCRSMatrix<B,A>::random) {
60
61 for (int i=0; i<size; i++)
63
65
66 for (int i=0; i<size; i++)
67 this->BCRSMatrix<B,A>::addindex(i, i);
68
70
71 }
72
74 BDMatrix (std::initializer_list<B> const &list)
75 : BDMatrix(list.size())
76 {
77 size_t i=0;
78 for (auto it = list.begin(); it != list.end(); ++it, ++i)
79 (*this)[i][i] = *it;
80 }
81
83 void setSize(size_type size)
84 {
85 this->BCRSMatrix<B,A>::setSize(size, // rows
86 size, // columns
87 size); // nonzeros
88
89 for (auto i : range(size))
91
93
94 for (auto i : range(size))
95 this->BCRSMatrix<B,A>::addindex(i, i);
96
98 }
99
101 BDMatrix& operator= (const BDMatrix& other) {
102 this->BCRSMatrix<B,A>::operator=(other);
103 return *this;
104 }
105
109 return *this;
110 }
111
117 template <class V>
118 void solve (V& x, const V& rhs) const {
119 for (size_type i=0; i<this->N(); i++)
120 {
121 auto&& xv = Impl::asVector(x[i]);
122 auto&& rhsv = Impl::asVector(rhs[i]);
123 Impl::asMatrix((*this)[i][i]).solve(xv,rhsv);
124 }
125 }
126
128 void invert() {
129 for (size_type i=0; i<this->N(); i++)
130 Impl::asMatrix((*this)[i][i]).invert();
131 }
132
133 private:
134
135 // ////////////////////////////////////////////////////////////////////////////
136 // The following methods from the base class should now actually be called
137 // ////////////////////////////////////////////////////////////////////////////
138
139 // createbegin and createend should be in there, too, but I can't get it to compile
140 // BCRSMatrix<B,A>::CreateIterator createbegin () {}
141 // BCRSMatrix<B,A>::CreateIterator createend () {}
142 void setrowsize (size_type i, size_type s) {}
143 void incrementrowsize (size_type i) {}
144 void endrowsizes () {}
145 void addindex (size_type row, size_type col) {}
146 void endindices () {}
147 };
150} // end namespace Dune
151
152#endif
Helper functions for determining the vector/matrix block level.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:464
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1147
void setrowsize(size_type i, size_type s)
Set number of indices in row i to s.
Definition: bcrsmatrix.hh:1115
@ 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
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:909
A block-diagonal matrix.
Definition: bdmatrix.hh:31
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bdmatrix.hh:37
A::size_type size_type
implement row_type with compressed vector
Definition: bdmatrix.hh:49
BDMatrix()
Default constructor.
Definition: bdmatrix.hh:56
BDMatrix(std::initializer_list< B > const &list)
Construct from a std::initializer_list.
Definition: bdmatrix.hh:74
B block_type
export the type representing the components
Definition: bdmatrix.hh:40
void solve(V &x, const V &rhs) const
Solve the system Ax=b in O(n) time.
Definition: bdmatrix.hh:118
A allocator_type
export the allocator type
Definition: bdmatrix.hh:43
BDMatrix & operator=(const BDMatrix &other)
assignment
Definition: bdmatrix.hh:101
static constexpr unsigned int blocklevel
increment block level counter
Definition: bdmatrix.hh:53
void setSize(size_type size)
Resize the matrix. Invalidates the content!
Definition: bdmatrix.hh:83
void invert()
Inverts the matrix.
Definition: bdmatrix.hh:128
Implementation of the BCRSMatrix class.
Dune namespace.
Definition: alignedallocator.hh:11
Utilities for reduction like operations on ranges.
Implements a scalar matrix view wrapper around an existing scalar.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)