Dune Core Modules (2.5.2)

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 
9 
15 namespace Dune {
25  template <class B, class A=std::allocator<B> >
26  class BDMatrix : public BCRSMatrix<B,A>
27  {
28  public:
29 
30  //===== type definitions and constants
31 
33  typedef typename B::field_type field_type;
34 
36  typedef B block_type;
37 
39  typedef A allocator_type;
40 
42  //typedef BCRSMatrix<B,A>::row_type row_type;
43 
45  typedef typename A::size_type size_type;
46 
48  enum {blocklevel = B::blocklevel+1};
49 
51  BDMatrix() : BCRSMatrix<B,A>() {}
52 
53  explicit BDMatrix(int size)
54  : BCRSMatrix<B,A>(size, size, BCRSMatrix<B,A>::random) {
55 
56  for (int i=0; i<size; i++)
57  this->BCRSMatrix<B,A>::setrowsize(i, 1);
58 
60 
61  for (int i=0; i<size; i++)
62  this->BCRSMatrix<B,A>::addindex(i, i);
63 
65 
66  }
67 
69  BDMatrix& operator= (const BDMatrix& other) {
70  this->BCRSMatrix<B,A>::operator=(other);
71  return *this;
72  }
73 
77  return *this;
78  }
79 
81  void invert() {
82  for (int i=0; i<this->N(); i++)
83  (*this)[i][i].invert();
84  }
85 
86  private:
87 
88  // ////////////////////////////////////////////////////////////////////////////
89  // The following methods from the base class should now actually be called
90  // ////////////////////////////////////////////////////////////////////////////
91 
92  // createbegin and createend should be in there, too, but I can't get it to compile
93  // BCRSMatrix<B,A>::CreateIterator createbegin () {}
94  // BCRSMatrix<B,A>::CreateIterator createend () {}
95  void setrowsize (size_type i, size_type s) {}
96  void incrementrowsize (size_type i) {}
97  void endrowsizes () {}
98  void addindex (size_type row, size_type col) {}
99  void endindices () {}
100  };
103 } // end namespace Dune
104 
105 #endif
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
BCRSMatrix & operator=(const BCRSMatrix &Mat)
assignment
Definition: bcrsmatrix.hh:878
void endrowsizes()
indicate that size of all rows is defined
Definition: bcrsmatrix.hh:1118
@ random
Build entries randomly.
Definition: bcrsmatrix.hh:489
void endindices()
indicate that all indices are defined, check consistency
Definition: bcrsmatrix.hh:1217
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1903
A block-diagonal matrix.
Definition: bdmatrix.hh:27
A::size_type size_type
implement row_type with compressed vector
Definition: bdmatrix.hh:45
BDMatrix()
Default constructor.
Definition: bdmatrix.hh:51
B block_type
export the type representing the components
Definition: bdmatrix.hh:36
A allocator_type
export the allocator type
Definition: bdmatrix.hh:39
BDMatrix & operator=(const BDMatrix &other)
assignment
Definition: bdmatrix.hh:69
B::field_type field_type
export the type representing the field
Definition: bdmatrix.hh:33
void invert()
Inverts the matrix.
Definition: bdmatrix.hh:81
Dune namespace.
Definition: alignment.hh:11
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 14, 22:30, 2024)