Newer
Older
#ifndef SUPPORT_PATCH_FACTORY_HH
#define SUPPORT_PATCH_FACTORY_HH
6
7
8
9
10
11
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
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
#include <dune/common/bitsetvector.hh>
#include <dune/common/fvector.hh>
#include <dune/fufem/referenceelementhelper.hh>
#include <dune/grid/common/mcmgmapper.hh>
#include <dune/faultnetworks/hierarchicleveliterator.hh>
template <class BasisType>
class SupportPatchFactory
{
protected:
typedef typename BasisType::GridView GV;
typedef typename BasisType::LocalFiniteElement LFE;
typedef typename GV::Grid GridType;
typedef typename GridType::ctype ctype;
static const int dim = GridType::dimension;
typedef typename GridType::template Codim<0>::Entity Element;
typedef typename GridType::template Codim<dim>::Entity Vertex;
typedef HierarchicLevelIterator<GridType> HierarchicLevelIteratorType;
typedef typename Dune::LevelMultipleCodimMultipleGeomTypeMapper<GridType, Dune::MCMGElementLayout > ElementMapper;
private:
const BasisType& coarseBasis_;
const BasisType& fineBasis_;
const int coarseLevel_;
const int fineLevel_;
const std::vector< std::vector <Element>>& vertexInElements_;
const int centerVertexIdx_;
const int patchDepth_;
const GridType& grid;
const GV& coarseGridView;
const GV& fineGridView;
std::vector<int> localToGlobal;
Dune::BitSetVector<1> boundaryDofs;
std::vector<Element> fineRegionElements;
ElementMapper coarseMapper;
public:
//setup
SupportPatchFactory(const BasisType& coarseBasis, const BasisType& fineBasis, const int coarseLevel, const int fineLevel, const std::vector< std::vector <Element>>& vertexInElements, const int centerVertexIdx, const int patchDepth) :
coarseBasis_(coarseBasis),
fineBasis_(fineBasis),
coarseLevel_(coarseLevel),
fineLevel_(fineLevel),
vertexInElements_(vertexInElements),
centerVertexIdx_(centerVertexIdx),
patchDepth_(patchDepth),
grid(coarseBasis_.getGridView().grid()),
coarseGridView(coarseBasis_.getGridView()),
fineGridView(fineBasis_.getGridView()),
coarseMapper(grid, coarseLevel_)
{
fineRegionElements.resize(0);
std::vector<Element> coarseElements;
Dune::BitSetVector<1> visited(coarseMapper.size());
Dune::BitSetVector<1> vertexVisited(vertexInElements_.size());
visited.unsetAll();
vertexVisited.unsetAll();
std::set<int> localDofs;
std::set<int> localBoundaryDofs;
addVertexSupport(centerVertexIdx_, coarseElements, visited, vertexVisited);
// construct coarse patch
for (int depth=1; depth<patchDepth_; depth++) {
std::vector<Element> coarseElementNeighbors;
for (size_t i=0; i<coarseElements.size(); ++i) {
const Element& elem = coarseElements[i];
const LFE& coarseLFE = coarseBasis_.getLocalFiniteElement(elem);
for (size_t j=0; j<coarseLFE.localBasis().size(); ++j) {
int dofIndex = coarseBasis_.indexInGridView(elem, j);
if (!vertexVisited[dofIndex][0]) {
addVertexSupport(dofIndex, coarseElementNeighbors, visited, vertexVisited);
}
}
}
coarseElements.insert(coarseElements.begin(), coarseElementNeighbors.begin(), coarseElementNeighbors.end());
}
// construct fine patch
for (size_t i=0; i<coarseElements.size(); ++i) {
const Element& elem = coarseElements[i];
addFinePatchElements(elem, visited, localDofs, localBoundaryDofs);
}
localToGlobal.resize(localDofs.size());
boundaryDofs.resize(localDofs.size());
boundaryDofs.unsetAll();
std::set<int>::iterator dofIt = localDofs.begin();
std::set<int>::iterator dofEndIt = localDofs.end();
size_t i=0;
for (; dofIt != dofEndIt; ++dofIt) {
int localDof = *dofIt;
localToGlobal[i] = localDof;
if (localBoundaryDofs.count(localDof)) {
boundaryDofs[i][0] = true;
}
i++;
}
}
std::vector<int>& getLocalToGlobal() {
return localToGlobal;
}
std::vector<Element>& getRegionElements() {
return fineRegionElements;
}
Dune::BitSetVector<1>& getBoundaryDofs() {
return boundaryDofs;
}
private:
void print(const Dune::BitSetVector<1>& x, std::string message){
std::cout << message << std::endl;
for (size_t i=0; i<x.size(); i++) {
std::cout << x[i][0] << " ";
}
std::cout << std::endl << std::endl;
}
void addVertexSupport(const int vertexIdx, std::vector<Element>& coarseElements, Dune::BitSetVector<1>& visited, Dune::BitSetVector<1>& vertexVisited) {
const std::vector<Element>& vertexInElement = vertexInElements_[vertexIdx];
for (size_t i=0; i<vertexInElement.size(); i++) {
const Element& elem = vertexInElement[i];
const int elemIdx = coarseMapper.index(elem);
if (!visited[elemIdx][0]) {
coarseElements.push_back(elem);
visited[elemIdx][0] = true;
}
}
vertexVisited[vertexIdx][0] = true;
}
bool containsInsideSubentity(const Element& elem, const typename GridType::LevelIntersection& intersection, int subEntity, int codim) {
return ReferenceElementHelper<double, dim>::subEntityContainsSubEntity(elem.type(), intersection.indexInInside(), 1, subEntity, codim);
}
void addFinePatchElements(const Element& coarseElem, const Dune::BitSetVector<1>& inCoarsePatch, std::set<int>& localDofs, std::set<int>& localBoundaryDofs) {
HierarchicLevelIteratorType endHierIt(grid, coarseElem, HierarchicLevelIteratorType::PositionFlag::end, fineLevel_);
for (HierarchicLevelIteratorType hierIt(grid, coarseElem, HierarchicLevelIteratorType::PositionFlag::begin, fineLevel_); hierIt!=endHierIt; ++hierIt) {
const Element& fineElem = *hierIt;
fineRegionElements.push_back(fineElem);
addLocalDofs(coarseElem, fineElem, inCoarsePatch, localDofs, localBoundaryDofs);
}
}
void addLocalDofs(const Element& coarseElem, const Element& fineElem, const Dune::BitSetVector<1>& inCoarsePatch, std::set<int>& localDofs, std::set<int>& localBoundaryDofs) {
const LFE& fineLFE = fineBasis_.getLocalFiniteElement(fineElem);
/*
const auto& fineGeometry = fineElem.geometry();
const auto& coarseGeometry = coarseElem.geometry();
const auto& ref = Dune::ReferenceElements<ctype, dim>::general(fineElem.type()); */
// insert local dofs
for (size_t i=0; i<fineLFE.localBasis().size(); ++i) {
int dofIndex = fineBasis_.index(fineElem, i);
localDofs.insert(dofIndex);
}
// search for parts of the global and inner boundary
typename GridType::LevelIntersectionIterator neighborIt = fineGridView.ibegin(fineElem);
typename GridType::LevelIntersectionIterator neighborItEnd = fineGridView.iend(fineElem);
for (; neighborIt!=neighborItEnd; ++neighborIt) {
const typename GridType::LevelIntersection& intersection = *neighborIt;
bool isInnerBoundary = false;
if (intersection.neighbor()) {
Element outsideFather = intersection.outside();
while (outsideFather.level()>coarseLevel_) {
outsideFather = outsideFather.father();
}
if (!inCoarsePatch[coarseMapper.index(outsideFather)][0]) {
isInnerBoundary = true;
}
}
if (intersection.boundary() or isInnerBoundary) {
typedef typename LFE::Traits::LocalCoefficientsType LocalCoefficients;
const LocalCoefficients* localCoefficients = &fineBasis_.getLocalFiniteElement(intersection.inside()).localCoefficients();
for (size_t i=0; i<localCoefficients->size(); i++) {
unsigned int entity = localCoefficients->localKey(i).subEntity();
unsigned int codim = localCoefficients->localKey(i).codim();
if (containsInsideSubentity(fineElem, intersection, entity, codim))
localBoundaryDofs.insert(fineBasis_.index(intersection.inside(), i));
}
}
}
}
};
#endif