AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
SimpleDataTransfer.hpp
1#pragma once
2
3#include <map>
4#include <optional>
5#include <vector>
6
7#include <amdis/DataTransfer.hpp>
8#include <amdis/functions/EntitySet.hpp>
9#include <amdis/gridfunctions/CoarsenedGridFunction.hpp>
10#include <amdis/gridfunctions/DiscreteFunction.hpp>
11
12#include <dune/functions/gridfunctions/finefunctiononcoarsegridview.hh>
13#include <dune/grid/common/capabilities.hh>
14
15namespace AMDiS {
16
23namespace tag {
24
29
30} // end namespace tag
31
32
33template <class Basis, class Container>
35{
36 using CoefficientVector = std::vector<typename Container::value_type>;
37 using EntityId = typename Basis::GridView::Grid::LocalIdSet::IdType;
38 using PersistentStorage = std::map<EntityId, CoefficientVector>;
39 using Domain = typename Basis::GridView::template Codim<0>::Entity::Geometry::LocalCoordinate;
40
41private:
42 PersistentStorage persistentStorage_;
43
44public:
46 {
47 static_assert(Dune::Capabilities::hasSingleGeometryType<typename Basis::GridView::Grid>::v);
48 }
49
51 void preAdapt(Basis const& basis, Container const& coeff, bool mightCoarsen)
52 {
53 CoefficientVector interpolationCoeff;
54 auto localView = basis.localView();
55 auto interpol = [&](auto const& lf, auto& localCoeff) {
56 localCoeff.resize(localView.size());
57 Traversal::forEachLeafNode(localView.tree(), [&](auto const& node, auto const& tp)
58 {
59 auto lf_at_tp = HierarchicNodeWrapper{tp,lf};
60 node.finiteElement().localInterpolation().interpolate(lf_at_tp, interpolationCoeff);
61 assert(node.size() == interpolationCoeff.size());
62 for (std::size_t i = 0; i < node.size(); ++i)
63 localCoeff[node.localIndex(i)] = interpolationCoeff[i];
64 });
65 };
66
67 DiscreteFunction coeffFct{coeff, basis};
68 auto c_gf = Dune::Functions::FineFunctionOnCoarseGridView{coeffFct, basis.gridView()};
69 auto c_lf = localFunction(c_gf);
70 auto const& idSet = basis.gridView().grid().localIdSet();
71 for (auto const& e : entitySet(basis))
72 {
73 localView.bind(e);
74
75 // store the coefficients on the current element
76 coeff.gather(localView, persistentStorage_[idSet.id(e)]);
77
78 // if the element might be coarsened, extract also interpolated coefficients on
79 // coarser elements
80 if (e.mightVanish()) {
81 auto father = e;
82 while (father.hasFather()) {
83 father = father.father();
84 auto [it,inserted] = persistentStorage_.try_emplace(idSet.id(father));
85 if (inserted) {
86 c_lf.bind(father);
87 interpol(c_lf, it->second);
88 }
89
90 if (!inserted || !father.mightVanish())
91 break;
92 }
93 }
94 }
95 }
96
97
99 void adapt(Basis const& basis, Container& coeff)
100 {
101 coeff.init(basis, false);
102 auto localView = basis.localView();
103
104 CoefficientVector interpolationCoeff;
105 auto interpol = [&](auto const& local, auto const& localCoeff) {
106 Traversal::forEachLeafNode(localView.tree(), [&](auto const& node, auto const& tp)
107 {
108 auto const& localFE = node.finiteElement();
109 auto const& localBasis = localFE.localBasis();
110 auto const& localInterpolation = localFE.localInterpolation();
111
112 using LocalBasis = TYPEOF(localBasis);
113 using RangeType = typename LocalBasis::Traits::RangeType;
114 std::vector<RangeType> shapeValues;
115 localInterpolation.interpolate([&](auto const& x) {
116 localBasis.evaluateFunction(local(x), shapeValues);
117 auto range = localCoeff[node.localIndex(0)] * shapeValues[0];
118 assert(node.size() == shapeValues.size());
119 for (std::size_t i = 1; i < shapeValues.size(); ++i)
120 range.axpy(localCoeff[node.localIndex(i)], shapeValues[i]);
121 return range;
122 }, interpolationCoeff);
123 coeff.scatter(localView, node, interpolationCoeff, Assigner::assign{});
124 });
125 };
126
127 auto const& idSet = basis.gridView().grid().localIdSet();
128 for (const auto& e : entitySet(basis))
129 {
130 localView.bind(e);
131 if (e.isNew()) {
132 // interpolate from coarser elements to the fine element
133 auto father = e;
134 std::function<Domain(Domain const&)> local = [](Domain const& x) { return x; };
135 [[maybe_unused]] bool found = false;
136 while (father.hasFather()) {
137 // store the coordinate transform
138 local = [local, geo=father.geometryInFather()](Domain const& x) {
139 return geo.global(local(x));
140 };
141 father = father.father();
142
143 // if coarse element is found, interpolate to fine element
144 auto it = persistentStorage_.find(idSet.id(father));
145 if (it != persistentStorage_.end()) {
146 interpol(local, it->second);
147 found = true;
148 break;
149 }
150 }
151 assert(found);
152 } else {
153 // just extract the stored coefficients
154 auto it = persistentStorage_.find(idSet.id(e));
155 test_exit(it != persistentStorage_.end(),
156 "Element should be in persistent storage but isn't.");
157 coeff.scatter(localView, it->second, Assigner::assign{});
158 }
159 }
160
161 coeff.finish();
162 }
163
165 void postAdapt(Container& coeff)
166 {
167 persistentStorage_.clear();
168 };
169};
170
172template <>
173struct DataTransferFactory<tag::simple_datatransfer>
174{
175 template <class Basis, class T, template <class> class Impl,
176 REQUIRES(Concepts::GlobalBasis<Basis>)>
177 static DataTransfer<Basis,VectorFacade<T,Impl>> create(Basis const&, VectorFacade<T,Impl> const&)
178 {
180 }
181};
182
184
185} // end namespace AMDiS
The base class for data transfer classes.
Definition DataTransfer.hpp:68
A mutable view on the subspace of a DOFVector,.
Definition DiscreteFunction.hpp:45
Definition SimpleDataTransfer.hpp:35
void preAdapt(Basis const &basis, Container const &coeff, bool mightCoarsen)
Create a persistent storage of the coefficients.
Definition SimpleDataTransfer.hpp:51
void postAdapt(Container &coeff)
Is called after the preAdapt step of the grid.
Definition SimpleDataTransfer.hpp:165
void adapt(Basis const &basis, Container &coeff)
Transfer the content of the coefficients to the new basis.
Definition SimpleDataTransfer.hpp:99
The basic container that stores a base vector and a corresponding basis.
Definition VectorFacade.hpp:49
Definition Assigner.hpp:8
Definition DataTransfer.hpp:133
Base tag type for all data transfer tags.
Definition DataTransfer.hpp:121
Tag that referrs to the SimpleDataTransfer.
Definition SimpleDataTransfer.hpp:28