SparsityPattern.hpp 3.19 KB
Newer Older
Praetorius, Simon's avatar
Praetorius, Simon committed
1
2
3
4
5
6
7
8
9
10
#pragma once

#include <memory>

#include <dune/istl/matrixindexset.hh>

#include <amdis/common/Index.hpp>

namespace AMDiS
{
11
  /// \brief A general sparsity pattern implementation using the full pattern of the
Praetorius, Simon's avatar
Praetorius, Simon committed
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
  /// basis by adding all local indices
  class SparsityPattern
  {
  public:
    /// Number of rows in the matrix
    std::size_t rows() const
    {
      return rows_;
    }

    /// Number of columns in the matrix
    std::size_t cols() const
    {
      return cols_;
    }

    /// Average number of non-zeros per row
    std::size_t avgRowSize() const
    {
      assert(rows_ == pattern_.rows());
      std::size_t sum = 0;
      for (std::size_t r = 0; r < rows_; ++r)
        sum += rowSize(r);
      return (sum + rows_ - 1) / rows_;
    }

    /// Number of non-zeros in row r
    std::size_t rowSize(std::size_t r) const
    {
      return pattern_.rowsize(r);
    }

    /// Total number of non-zeros
    std::size_t nnz() const
    {
      return pattern_.size();
    }

    /// Initialize a matrix with the non-zero structure
    template <class Matrix>
    void applyTo(Matrix& matrix) const
    {
      pattern_.exportIdx(matrix);
    }

    // Update pattern when basis is updated. This method is called if rowBasis == colBasis.
58
59
    template <class Basis>
    void init(Basis const& basis)
Praetorius, Simon's avatar
Praetorius, Simon committed
60
    {
61
      rows_ = basis.dimension();
Praetorius, Simon's avatar
Praetorius, Simon committed
62
      cols_ = rows_;
63
      pattern_.resize(0, 0); // clear the old pattern
Praetorius, Simon's avatar
Praetorius, Simon committed
64
65
      pattern_.resize(rows_, cols_);

66
67
      auto localView = basis.localView();
      for (const auto& element : elements(basis.gridView())) {
Praetorius, Simon's avatar
Praetorius, Simon committed
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
        localView.bind(element);

        for (std::size_t i = 0, size = localView.size(); i < size; ++i) {
          // The global index of the i-th vertex of the element
          auto row = localView.index(i);
          for (std::size_t j = 0; j < size; ++j) {
            // The global index of the j-th vertex of the element
            auto col = localView.index(j);
            pattern_.add(row, col);
          }
        }
      }
    }

    // Update pattern when basis is updated. This method is called if rowBasis != colBasis.
83
84
    template <class RowBasis, class ColBasis>
    void init(RowBasis const& rowBasis, ColBasis const& colBasis)
Praetorius, Simon's avatar
Praetorius, Simon committed
85
    {
86
87
88
89
90
91
      if (uintptr_t(&rowBasis) == uintptr_t(&colBasis))
        return init(rowBasis);

      rows_ = rowBasis.dimension();
      cols_ = colBasis.dimension();
      pattern_.resize(0, 0); // clear the old pattern
Praetorius, Simon's avatar
Praetorius, Simon committed
92
93
      pattern_.resize(rows_, cols_);

94
95
96
      auto rowLocalView = rowBasis.localView();
      auto colLocalView = colBasis.localView();
      for (const auto& element : elements(rowBasis.gridView())) {
Praetorius, Simon's avatar
Praetorius, Simon committed
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
        rowLocalView.bind(element);
        colLocalView.bind(element);

        for (std::size_t i = 0; i < rowLocalView.size(); ++i) {
          // The global index of the i-th vertex of the element
          auto row = rowLocalView.index(i);
          for (std::size_t j = 0; j < colLocalView.size(); ++j) {
            // The global index of the j-th vertex of the element
            auto col = colLocalView.index(j);
            pattern_.add(row, col);
          }
        }
      }
    }

  private:
    std::size_t rows_;
    std::size_t cols_;
    Dune::MatrixIndexSet pattern_;
  };

} // end namespace AMDiS