Dune Core Modules (2.9.0)

bdmatrix.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_BDMATRIX_HH
6 #define DUNE_ISTL_BDMATRIX_HH
7 
8 #include <memory>
9 
12 
13 #include <dune/istl/bcrsmatrix.hh>
14 #include <dune/istl/blocklevel.hh>
15 
21 namespace Dune {
31  template <class B, class A=std::allocator<B> >
32  class BDMatrix : public BCRSMatrix<B,A>
33  {
34  public:
35 
36  //===== type definitions and constants
37 
39  using field_type = typename Imp::BlockTraits<B>::field_type;
40 
42  typedef B block_type;
43 
45  typedef A allocator_type;
46 
48  //typedef BCRSMatrix<B,A>::row_type row_type;
49 
51  typedef typename A::size_type size_type;
52 
54  [[deprecated("Use free function blockLevel(). Will be removed after 2.8.")]]
55  static constexpr unsigned int blocklevel = blockLevel<B>()+1;
56 
58  BDMatrix() : BCRSMatrix<B,A>() {}
59 
60  explicit BDMatrix(int size)
61  : BCRSMatrix<B,A>(size, size, BCRSMatrix<B,A>::random) {
62 
63  for (int i=0; i<size; i++)
64  this->BCRSMatrix<B,A>::setrowsize(i, 1);
65 
67 
68  for (int i=0; i<size; i++)
69  this->BCRSMatrix<B,A>::addindex(i, i);
70 
72 
73  }
74 
76  BDMatrix (std::initializer_list<B> const &list)
77  : BDMatrix(list.size())
78  {
79  size_t i=0;
80  for (auto it = list.begin(); it != list.end(); ++it, ++i)
81  (*this)[i][i] = *it;
82  }
83 
85  void setSize(size_type size)
86  {
87  this->BCRSMatrix<B,A>::setSize(size, // rows
88  size, // columns
89  size); // nonzeros
90 
91  for (auto i : range(size))
92  this->BCRSMatrix<B,A>::setrowsize(i, 1);
93 
95 
96  for (auto i : range(size))
97  this->BCRSMatrix<B,A>::addindex(i, i);
98 
100  }
101 
103  BDMatrix& operator= (const BDMatrix& other) {
104  this->BCRSMatrix<B,A>::operator=(other);
105  return *this;
106  }
107 
111  return *this;
112  }
113 
119  template <class V>
120  void solve (V& x, const V& rhs) const {
121  for (size_type i=0; i<this->N(); i++)
122  {
123  auto&& xv = Impl::asVector(x[i]);
124  auto&& rhsv = Impl::asVector(rhs[i]);
125  Impl::asMatrix((*this)[i][i]).solve(xv,rhsv);
126  }
127  }
128 
130  void invert() {
131  for (size_type i=0; i<this->N(); i++)
132  Impl::asMatrix((*this)[i][i]).invert();
133  }
134 
135  private:
136 
137  // ////////////////////////////////////////////////////////////////////////////
138  // The following methods from the base class should now actually be called
139  // ////////////////////////////////////////////////////////////////////////////
140 
141  // createbegin and createend should be in there, too, but I can't get it to compile
142  // BCRSMatrix<B,A>::CreateIterator createbegin () {}
143  // BCRSMatrix<B,A>::CreateIterator createend () {}
144  void setrowsize (size_type i, size_type s) {}
145  void incrementrowsize (size_type i) {}
146  void endrowsizes () {}
147  void addindex (size_type row, size_type col) {}
148  void endindices () {}
149  };
150 
151  template<typename B, typename A>
152  struct FieldTraits< BDMatrix<B, A> >
153  {
154  using field_type = typename BDMatrix<B, A>::field_type;
155  using real_type = typename FieldTraits<field_type>::real_type;
156  };
159 } // end namespace Dune
160 
161 #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
void setrowsize(size_type i, size_type s)
Set number of indices in row i to s.
Definition: bcrsmatrix.hh:1117
@ 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-diagonal matrix.
Definition: bdmatrix.hh:33
typename Imp::BlockTraits< B >::field_type field_type
export the type representing the field
Definition: bdmatrix.hh:39
A::size_type size_type
implement row_type with compressed vector
Definition: bdmatrix.hh:51
BDMatrix()
Default constructor.
Definition: bdmatrix.hh:58
BDMatrix(std::initializer_list< B > const &list)
Construct from a std::initializer_list.
Definition: bdmatrix.hh:76
B block_type
export the type representing the components
Definition: bdmatrix.hh:42
void solve(V &x, const V &rhs) const
Solve the system Ax=b in O(n) time.
Definition: bdmatrix.hh:120
A allocator_type
export the allocator type
Definition: bdmatrix.hh:45
BDMatrix & operator=(const BDMatrix &other)
assignment
Definition: bdmatrix.hh:103
static constexpr unsigned int blocklevel
increment block level counter
Definition: bdmatrix.hh:55
void setSize(size_type size)
Resize the matrix. Invalidates the content!
Definition: bdmatrix.hh:85
void invert()
Inverts the matrix.
Definition: bdmatrix.hh:130
Dune namespace.
Definition: alignedallocator.hh:13
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.80.0 (May 4, 22:30, 2024)