3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
partialreassembler.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_PARTIAL_REASSEMBLER_HH
25#define DUMUX_PARTIAL_REASSEMBLER_HH
26
27#include <algorithm>
28#include <vector>
29
30#include <dune/grid/common/gridenums.hh>
31#include <dune/istl/multitypeblockvector.hh>
32
36
37#include "entitycolor.hh"
38
39namespace Dumux {
40
42{
43public:
44 template<typename... Args>
46 { DUNE_THROW(Dune::InvalidStateException, "DefaultPartialReassembler should be never constructed!"); }
47
48 template<typename... Args>
49 void report(Args&&... args) {}
50
51 template<typename... Args>
52 void resetJacobian(Args&&... args) const {}
53
54 template<typename... Args>
55 void computeColors(Args&&... args) {}
56
57 template<typename... Args>
58 void resetColors(Args&&... args) {}
59
60 EntityColor dofColor(size_t idx) const
61 { return EntityColor::red; }
62
63 EntityColor elementColor(size_t idx) const
64 { return EntityColor::red; }
65
66 EntityColor vertexColor(size_t idx) const
67 { return EntityColor::red; }
68};
69
71template<class Assembler, DiscretizationMethod discMethod>
73{
74public:
75 PartialReassemblerEngine(const Assembler&)
76 { DUNE_THROW(Dune::NotImplemented, "PartialReassembler for this discretization method!"); }
77
78 EntityColor elementColor(size_t idx) const
79 { return EntityColor::red; }
80
81 EntityColor dofColor(size_t idx) const
82 { return EntityColor::red; }
83
84 template<typename... Args>
85 std::size_t computeColors(Args&&... args) { return 0; }
86
87 template<typename... Args>
88 void resetJacobian(Args&&... args) const {}
89
90 template<typename... Args>
91 void resetColors(Args&&... args) {}
92};
93
98template<class Assembler>
100{
101 using Scalar = typename Assembler::Scalar;
102 using GridGeometry = typename Assembler::GridGeometry;
103 using JacobianMatrix = typename Assembler::JacobianMatrix;
104 using VertexMapper = typename GridGeometry::VertexMapper;
105 static constexpr int dim = GridGeometry::GridView::dimension;
106
107public:
108 PartialReassemblerEngine(const Assembler& assembler)
109 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
110 , vertexColor_(assembler.gridGeometry().vertexMapper().size(), EntityColor::red)
111 {}
112
113 // returns number of green elements
114 std::size_t computeColors(const Assembler& assembler,
115 const std::vector<Scalar>& distanceFromLastLinearization,
116 Scalar threshold)
117 {
118 const auto& gridGeometry = assembler.gridGeometry();
119 const auto& gridView = gridGeometry.gridView();
120 const auto& elementMapper = gridGeometry.elementMapper();
121 const auto& vertexMapper = gridGeometry.vertexMapper();
122
123 // set all vertices to green
124 vertexColor_.assign(vertexColor_.size(), EntityColor::green);
125
126 // mark the red vertices
127 for (unsigned int i = 0; i < vertexColor_.size(); ++i)
128 {
129 using std::max;
130 if (distanceFromLastLinearization[i] > threshold)
131 // mark vertex as red if discrepancy is larger than
132 // the relative tolerance
133 vertexColor_[i] = EntityColor::red;
134 }
135
136 // Mark all red elements
137 for (const auto& element : elements(gridView))
138 {
139 // find out whether the current element features a red vertex
140 bool isRed = false;
141
142 int numVertices = element.subEntities(dim);
143
144 for (int i = 0; i < numVertices; ++i) {
145 int globalI = vertexMapper.subIndex(element, i, dim);
146
147 if (vertexColor_[globalI] == EntityColor::red) {
148 isRed = true;
149 break;
150 }
151 }
152
153 int eIdx = elementMapper.index(element);
154 // if a vertex is red, the element color is also red, otherwise green
155 if (isRed)
156 elementColor_[eIdx] = EntityColor::red;
157 else
158 elementColor_[eIdx] = EntityColor::green;
159 }
160
161 // mark orange vertices
162 for (const auto& element : elements(gridView))
163 {
164 int eIdx = elementMapper.index(element);
165
166 // only red elements tint vertices yellow
167 if (elementColor_[eIdx] == EntityColor::red)
168 {
169 int numVertices = element.subEntities(dim);
170
171 for (int i = 0; i < numVertices; ++i) {
172 int globalI = vertexMapper.subIndex(element, i, dim);
173
174 // red vertices don't become orange
175 if (vertexColor_[globalI] != EntityColor::red)
176 vertexColor_[globalI] = EntityColor::orange;
177 }
178 }
179 }
180
181 // at this point we communicate the yellow vertices to the
182 // neighboring processes because a neigbor process may not see
183 // the red vertex for yellow border vertices
185 minHandle(vertexMapper, vertexColor_);
186 gridView.communicate(minHandle,
187 Dune::InteriorBorder_InteriorBorder_Interface,
188 Dune::ForwardCommunication);
189
190 // mark yellow elements
191 for (const auto& element : elements(gridView))
192 {
193 int eIdx = elementMapper.index(element);
194
195 // only treat non-red elements
196 if (elementColor_[eIdx] != EntityColor::red)
197 {
198 // check whether the element features a orange vertex
199 bool isOrange = false;
200 int numVertices = element.subEntities(dim);
201
202 for (int i = 0; i < numVertices; ++i) {
203 int globalI = vertexMapper.subIndex(element, i, dim);
204
205 if (vertexColor_[globalI] == EntityColor::orange) {
206 isOrange = true;
207 break;
208 }
209 }
210
211 if (isOrange)
212 elementColor_[eIdx] = EntityColor::yellow;
213 }
214 }
215
216 // change orange vertices to yellow ones if it has at least
217 // one green element as a neighbor
218 for (const auto& element : elements(gridView))
219 {
220 int eIdx = elementMapper.index(element);
221
222 // only green elements are considered
223 if (elementColor_[eIdx] == EntityColor::green)
224 {
225 int numVertices = element.subEntities(dim);
226
227 for (int i = 0; i < numVertices; ++i) {
228 int globalI = vertexMapper.subIndex(element, i, dim);
229
230 // if a vertex is orange, recolor it to yellow
231 if (vertexColor_[globalI] == EntityColor::orange)
232 vertexColor_[globalI] = EntityColor::yellow;
233 }
234 }
235 }
236
237 // demote the border orange vertices
239 maxHandle(vertexMapper, vertexColor_);
240 gridView.communicate(maxHandle,
241 Dune::InteriorBorder_InteriorBorder_Interface,
242 Dune::ForwardCommunication);
243
244 // promote the remaining orange vertices to red
245 for (unsigned int i=0; i < vertexColor_.size(); ++i) {
246 // if a vertex is green or yellow don't do anything
247 if (vertexColor_[i] == EntityColor::green || vertexColor_[i] == EntityColor::yellow)
248 continue;
249
250 // set the vertex to red
251 vertexColor_[i] = EntityColor::red;
252 }
253
254 // count green elements
255 return std::count_if(elementColor_.begin(), elementColor_.end(),
256 [](EntityColor c){ return c == EntityColor::green; });
257 }
258
259 void resetJacobian(Assembler& assembler) const
260 {
261 auto& jacobian = assembler.jacobian();
262
263 // loop over all dofs
264 for (unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
265 {
266 // reset all entries corrosponding to a non-green vertex
267 if (vertexColor_[rowIdx] != EntityColor::green)
268 {
269 // set all matrix entries in the row to 0
270 auto colIt = jacobian[rowIdx].begin();
271 const auto& colEndIt = jacobian[rowIdx].end();
272 for (; colIt != colEndIt; ++colIt) {
273 *colIt = 0.0;
274 }
275 }
276 }
277 }
278
280 {
281 elementColor_.assign(elementColor_.size(), EntityColor::red);
282 vertexColor_.assign(vertexColor_.size(), EntityColor::red);
283 }
284
285 EntityColor elementColor(size_t idx) const
286 { return elementColor_[idx]; }
287
288 EntityColor vertexColor(size_t idx) const
289 { return vertexColor_[idx]; }
290
291 EntityColor dofColor(size_t idx) const
292 { return vertexColor_[idx]; }
293
294private:
296 std::vector<EntityColor> elementColor_;
297 std::vector<EntityColor> vertexColor_;
298};
299
304template<class Assembler>
306{
307 using Scalar = typename Assembler::Scalar;
308 using GridGeometry = typename Assembler::GridGeometry;
309 using JacobianMatrix = typename Assembler::JacobianMatrix;
310
311public:
312 PartialReassemblerEngine(const Assembler& assembler)
313 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
314 {}
315
316 // returns number of green elements
317 std::size_t computeColors(const Assembler& assembler,
318 const std::vector<Scalar>& distanceFromLastLinearization,
319 Scalar threshold)
320 {
321 const auto& gridGeometry = assembler.gridGeometry();
322 const auto& gridView = gridGeometry.gridView();
323 const auto& elementMapper = gridGeometry.elementMapper();
324
325 // mark the red elements
326 for (const auto& element : elements(gridView))
327 {
328 int eIdx = elementMapper.index(element);
329
330 if (distanceFromLastLinearization[eIdx] > threshold)
331 {
332 // mark element as red if discrepancy is larger than
333 // the relative tolerance
334 elementColor_[eIdx] = EntityColor::red;
335 }
336 else
337 {
338 elementColor_[eIdx] = EntityColor::green;
339 }
340 }
341
342 // mark the neighbors also red
343 const auto& connectivityMap = gridGeometry.connectivityMap();
344 for (unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
345 {
346 if (elementColor_[eIdx] == EntityColor::red)
347 continue; // element is red already!
348
349 if (distanceFromLastLinearization[eIdx] > threshold)
350 {
351 for (const auto& connectedDof : connectivityMap[eIdx])
352 elementColor_[connectedDof.globalJ] = EntityColor::red;
353 }
354 }
355
356 // count green elements
357 return std::count_if(elementColor_.begin(), elementColor_.end(),
358 [](EntityColor c){return c == EntityColor::green;});
359
360 }
361
362 void resetJacobian(Assembler& assembler) const
363 {
364 auto& jacobian = assembler.jacobian();
365 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
366
367 // loop over all dofs
368 for (unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
369 {
370 // reset all entries corresponding to a non-green element
371 if (elementColor_[colIdx] != EntityColor::green)
372 {
373 // set all matrix entries in the column to 0
374 jacobian[colIdx][colIdx] = 0;
375 for (const auto& dataJ : connectivityMap[colIdx])
376 jacobian[dataJ.globalJ][colIdx] = 0;
377 }
378 }
379 }
380
382 {
383 elementColor_.assign(elementColor_.size(), EntityColor::red);
384 }
385
386 EntityColor elementColor(size_t idx) const
387 { return elementColor_[idx]; }
388
389 EntityColor dofColor(size_t idx) const
390 { return elementColor_[idx]; }
391
392private:
394 std::vector<EntityColor> elementColor_;
395};
396
401template<class Assembler>
403: public PartialReassemblerEngine<Assembler, DiscretizationMethod::cctpfa>
404{
406public:
407 using ParentType::ParentType;
408};
409
412{
413 template<class Engine>
414 auto operator()(Engine&& e) -> decltype(e.vertexColor(0)) {}
415};
416
422template<class Assembler>
424{
425 using Scalar = typename Assembler::Scalar;
426 using GridGeometry = typename Assembler::GridGeometry;
427 using JacobianMatrix = typename Assembler::JacobianMatrix;
428 using VertexMapper = typename GridGeometry::VertexMapper;
429
430 static constexpr DiscretizationMethod discMethod = GridGeometry::discMethod;
432
433public:
434
439 PartialReassembler(const Assembler& assembler)
440 : engine_(assembler)
441 , greenElems_(0)
442 {
443 const auto& gridGeometry = assembler.gridGeometry();
444 totalElems_ = gridGeometry.elementMapper().size();
445 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
446 }
447
463 void computeColors(const Assembler& assembler,
464 const std::vector<Scalar>& distanceFromLastLinearization,
465 Scalar threshold)
466 {
467 greenElems_ = engine_.computeColors(assembler, distanceFromLastLinearization, threshold);
468 }
469
471 {
472 engine_.resetColors();
473 }
474
475 void resetJacobian(Assembler& assembler) const
476 {
477 engine_.resetJacobian(assembler);
478 }
479
483 template <class Communication>
484 void report(const Communication& comm, std::ostream& outStream)
485 {
486 if (comm.size() > 1)
487 greenElems_ = comm.sum(greenElems_);
488
489 auto reassembledElems = totalElems_ - greenElems_;
490 auto width = std::to_string(totalElems_).size();
491 outStream << ", reassembled " << std::setw(width)
492 << reassembledElems << " (" << std::setw(3)
493 << 100*reassembledElems/totalElems_ << "%) elements";
494 }
495
496 EntityColor elementColor(size_t idx) const
497 { return engine_.elementColor(idx); }
498
499 EntityColor dofColor(size_t idx) const
500 { return engine_.dofColor(idx); }
501
502 template<bool enable = decltype(isValid(hasVertexColor()).template check<Engine>())::value,
503 typename std::enable_if_t<enable, int> = 0>
504 EntityColor vertexColor(size_t idx) const
505 { return engine_.vertexColor(idx); }
506
507private:
508 Engine engine_;
509 size_t totalElems_;
510 size_t greenElems_;
511};
512
513} // namespace Dumux
514
515#endif
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
A helper function for class member function introspection.
The available discretization methods in Dumux.
Contains a class to exchange entries of a vector.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
EntityColor
The colors of elements and vertices required for partial Jacobian reassembly.
Definition: entitycolor.hh:35
@ green
does not need to be reassembled
@ red
distance from last linearization is above the tolerance
@ yellow
neighboring entity is red
Definition: adapt.hh:29
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
Definition: partialreassembler.hh:42
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:60
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:52
void computeColors(Args &&... args)
Definition: partialreassembler.hh:55
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:66
void report(Args &&... args)
Definition: partialreassembler.hh:49
DefaultPartialReassembler(Args &&... args)
Definition: partialreassembler.hh:45
void resetColors(Args &&... args)
Definition: partialreassembler.hh:58
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:63
the partial reassembler engine specialized for discretization methods
Definition: partialreassembler.hh:73
std::size_t computeColors(Args &&... args)
Definition: partialreassembler.hh:85
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:78
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:81
void resetColors(Args &&... args)
Definition: partialreassembler.hh:91
PartialReassemblerEngine(const Assembler &)
Definition: partialreassembler.hh:75
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:88
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:291
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:108
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:288
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:285
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:114
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:259
The partial reassembler engine specialized for the cellcentered TPFA method.
Definition: partialreassembler.hh:306
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:362
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:386
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:312
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:317
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:389
helper struct to determine whether the an engine class has vertex colors
Definition: partialreassembler.hh:412
auto operator()(Engine &&e) -> decltype(e.vertexColor(0))
Definition: partialreassembler.hh:414
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:424
PartialReassembler(const Assembler &assembler)
constructor
Definition: partialreassembler.hh:439
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:504
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:475
void resetColors()
Definition: partialreassembler.hh:470
void computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Determine the colors of entities for partial reassembly.
Definition: partialreassembler.hh:463
void report(const Communication &comm, std::ostream &outStream)
called by the assembler after successful assembly
Definition: partialreassembler.hh:484
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:499
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:496
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:79