DOFVector.hpp 6.36 KB
Newer Older
1
2
3
4
5
#pragma once

#include <memory>
#include <utility>

6
#include <dune/common/concept.hh>
7
8
#include <dune/common/shared_ptr.hh>

9
10
#include <dune/functions/functionspacebases/concepts.hh>

11
12
#include <amdis/DataTransfer.hpp>
#include <amdis/LinearAlgebra.hpp>
13
#include <amdis/Observer.hpp>
14
15
#include <amdis/common/Concepts.hpp>
#include <amdis/common/TypeTraits.hpp>
16
#include <amdis/functions/ParallelGlobalBasis.hpp>
17
18
19
20
#include <amdis/typetree/TreePath.hpp>

namespace AMDiS
{
21
22
23
24
25
  // Forward declarations
  template <class PB>
  class ParallelGlobalBasis;


26
27
28
29
  /// \brief The basic container that stores a base vector and a corresponding basis
  /**
   * \tparam GB  Basis of the vector
   * \tparam T   Type of the coefficients
Praetorius, Simon's avatar
Praetorius, Simon committed
30
   * \tparam TraitsType  Collection of parameter for the linear-algebra backend
31
   **/
Praetorius, Simon's avatar
Praetorius, Simon committed
32
  template <class GB, class T = double, class TraitsType = BackendTraits<GB>>
33
  class DOFVector
Praetorius, Simon's avatar
Praetorius, Simon committed
34
      : public VectorFacade<T, TraitsType::template VectorImpl>
Praetorius, Simon's avatar
Praetorius, Simon committed
35
36
37
      , private Observer<event::preAdapt>
      , private Observer<event::adapt>
      , private Observer<event::postAdapt>
38
39
40
41
  {
    using Self = DOFVector;

  public:
42
    /// A global basis associated to the coefficients
43
44
    using GlobalBasis = GB;

45
46
47
    /// The internal coefficient vector storage
    using Coefficients = VectorFacade<T, TraitsType::template VectorImpl>;

48
49
50
51
    /// The index/size - type
    using size_type  = typename GlobalBasis::size_type;

    /// The type of the elements of the DOFVector
Praetorius, Simon's avatar
Praetorius, Simon committed
52
    using value_type = T;
53

54
55
    /// Wrapper for the implementation of the transfer of data during grid adaption
    using DataTransfer = DataTransferWrapper<Self>;
56

Praetorius, Simon's avatar
Praetorius, Simon committed
57
58
59
    /// Collection of parameter for the linear-algebra backend
    using Traits = TraitsType;

60
  public:
61
    /// (1) Constructor. Stores the shared_ptr of the basis and creates a new DataTransfer.
Praetorius, Simon's avatar
Praetorius, Simon committed
62
    DOFVector(std::shared_ptr<GB> const& basis,
63
              DataTransferOperation op = DataTransferOperation::INTERPOLATE)
Praetorius, Simon's avatar
Praetorius, Simon committed
64
      : Coefficients(*basis)
Praetorius, Simon's avatar
Praetorius, Simon committed
65
      , Observer<event::preAdapt>(basis->gridView().grid())
66
      , Observer<event::adapt>(*basis)
Praetorius, Simon's avatar
Praetorius, Simon committed
67
      , Observer<event::postAdapt>(basis->gridView().grid())
68
69
      , dataTransfer_(op, basis)
      , basis_(basis)
70
71
    {}

72
73
    /// (2) Constructor. Forwards to (1) by wrapping reference into non-destroying
    /// shared_ptr, see \ref Dune::wrap_or_move.
74
75
76
77
78
79
    template <class GB_,
      REQUIRES(Concepts::Similar<GB_,GB>)>
    DOFVector(GB_&& basis, DataTransferOperation op = DataTransferOperation::INTERPOLATE)
      : DOFVector(Dune::wrap_or_move(FWD(basis)), op)
    {}

80
81
82
83
84
    /// (3) Constructor. Forwards to (1) by creating a new basis from a dune-functions pre-basis factory.
    template <class GV, class PBF,
      REQUIRES(Concepts::PreBasisFactory<PBF, GV, MultiIndex_t<PBF>>)>
    DOFVector(GV const& gridView, PBF const& preBasisFactory, DataTransferOperation op = DataTransferOperation::INTERPOLATE)
      : DOFVector(std::make_shared<GB>(gridView, preBasisFactory), op)
85
86
    {}

87
    /// Return the global basis
88
    std::shared_ptr<GlobalBasis const> const& basis() const
89
    {
90
      return basis_;
91
92
    }

93
    Coefficients const& coefficients() const
Praetorius, Simon's avatar
Praetorius, Simon committed
94
95
96
97
    {
      return static_cast<Coefficients const&>(*this);
    }

98
    Coefficients& coefficients()
Praetorius, Simon's avatar
Praetorius, Simon committed
99
100
101
102
    {
      return static_cast<Coefficients&>(*this);
    }

103
104
105
106
107
108
    /// Write DOFVector to file
    void backup(std::string const& filename);

    /// Read backup data from file
    void restore(std::string const& filename);

109
110
111
112
113
114
115
116
117
118
119
120
121
    void resize()
    {
      Coefficients::resize(sizeInfo(*basis_));
    }

    using Coefficients::resize;

    void resizeZero()
    {
      Coefficients::resizeZero(sizeInfo(*basis_));
    }

    using Coefficients::resizeZero;
122
123

    /// Return the associated DataTransfer object
124
    DataTransfer const& dataTransfer() const
125
126
127
128
129
    {
      return dataTransfer_;
    }

    /// Return the associated DataTransfer object
130
    DataTransfer& dataTransfer()
131
132
133
134
135
    {
      return dataTransfer_;
    }

    /// Assign the DataTransfer object
136
    void setDataTransfer(DataTransfer const& dataTransfer)
137
138
139
140
    {
      dataTransfer_ = dataTransfer;
    }

141
142
143
144
    void setDataTransfer(DataTransfer&& dataTransfer)
    {
      dataTransfer_ = std::move(dataTransfer);
    }
145

146
147
    /// Create a new DataTransfer object based on the operation type
    void setDataTransfer(DataTransferOperation op)
148
    {
149
      setDataTransfer(DataTransfer(op, basis_));
150
151
    }

Praetorius, Simon's avatar
Praetorius, Simon committed
152
153
  protected:

154
155
    /// Override of Observer update(event::preAdapt) method. Redirects to preAdapt method of the
    /// \ref DataTransfer object.
Praetorius, Simon's avatar
Praetorius, Simon committed
156
    void updateImpl(event::preAdapt e) override
157
    {
Praetorius, Simon's avatar
Praetorius, Simon committed
158
      dataTransfer_.preAdapt(*this, e.value);
159
160
    }

161
162
    /// Override of Observer update(event::adapt) method. Redirects to adapt method of the
    /// \ref DataTransfer object.
Praetorius, Simon's avatar
Praetorius, Simon committed
163
    void updateImpl(event::adapt e) override
164
    {
Praetorius, Simon's avatar
Praetorius, Simon committed
165
      assert(e.value);
166
      resize();
167
      dataTransfer_.adapt(*this);
168
169
    }

170
171
    /// Override of Observer update(event::postAdapt) method. Redirects to postAdapt method of the
    /// \ref DataTransfer object.
Praetorius, Simon's avatar
Praetorius, Simon committed
172
    void updateImpl(event::postAdapt) override
173
    {
174
      dataTransfer_.postAdapt(*this);
175
176
177
178
179
    }

  private:
    /// Data interpolation when the grid changes, set by default
    /// to \ref DataTransferOperation::INTERPOLATE.
180
181
    DataTransfer dataTransfer_;

Praetorius, Simon's avatar
Praetorius, Simon committed
182
    std::shared_ptr<GlobalBasis const> basis_;
183
184
  };

185
186

  // deduction guides
187
188
  template <class GB>
  DOFVector(GB&& basis, DataTransferOperation op = DataTransferOperation::INTERPOLATE)
189
    -> DOFVector<Underlying_t<GB>>;
190

191
192
193
194
  template <class GV, class PBF>
  DOFVector(GV const& gridView, PBF const& pbf, DataTransferOperation op = DataTransferOperation::INTERPOLATE)
    -> DOFVector<decltype(ParallelGlobalBasis{gridView,pbf})>;

195
196
197
198
199
200
201
202
203
204
205

  /// \brief Create a DOFVector from a basis.
  /**
   * This generator function accepts the basis as reference, temporary, or
   * shared_ptr. Internally the reference is wrapped into a non-destroying
   * shared_ptr and the temporary is moved into a new shared_ptr.
   *
   * The DataTransferOperation controls what is done during grid changes with the
   * DOFVector. The default is interpolation of the data to the new grid. See
   * \ref DataTransferOperation for more options.
   **/
206
  template <class ValueType = double, class GB>
207
  DOFVector<Underlying_t<GB>, ValueType>
208
  makeDOFVector(GB&& basis, DataTransferOperation op = DataTransferOperation::INTERPOLATE)
209
210
211
212
213
214
215
  {
    return {FWD(basis), op};
  }

} // end namespace AMDiS

#include <amdis/DOFVector.inc.hpp>