3.3.0
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
33#include <dumux/io/format.hh>
37
38#include "entitycolor.hh"
39
40namespace Dumux {
41
43{
44public:
45 template<typename... Args>
47 { DUNE_THROW(Dune::InvalidStateException, "DefaultPartialReassembler should be never constructed!"); }
48
49 template<typename... Args>
50 void report(Args&&... args) {}
51
52 template<typename... Args>
53 void resetJacobian(Args&&... args) const {}
54
55 template<typename... Args>
56 void computeColors(Args&&... args) {}
57
58 template<typename... Args>
59 void resetColors(Args&&... args) {}
60
61 EntityColor dofColor(size_t idx) const
62 { return EntityColor::red; }
63
64 EntityColor elementColor(size_t idx) const
65 { return EntityColor::red; }
66
67 EntityColor vertexColor(size_t idx) const
68 { return EntityColor::red; }
69};
70
72template<class Assembler, DiscretizationMethod discMethod>
74{
75public:
76 PartialReassemblerEngine(const Assembler&)
77 { DUNE_THROW(Dune::NotImplemented, "PartialReassembler for this discretization method!"); }
78
79 EntityColor elementColor(size_t idx) const
80 { return EntityColor::red; }
81
82 EntityColor dofColor(size_t idx) const
83 { return EntityColor::red; }
84
85 template<typename... Args>
86 std::size_t computeColors(Args&&... args) { return 0; }
87
88 template<typename... Args>
89 void resetJacobian(Args&&... args) const {}
90
91 template<typename... Args>
92 void resetColors(Args&&... args) {}
93};
94
99template<class Assembler>
101{
102 using Scalar = typename Assembler::Scalar;
103 using GridGeometry = typename Assembler::GridGeometry;
104 using JacobianMatrix = typename Assembler::JacobianMatrix;
105 using VertexMapper = typename GridGeometry::VertexMapper;
106 static constexpr int dim = GridGeometry::GridView::dimension;
107
108public:
109 PartialReassemblerEngine(const Assembler& assembler)
110 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
111 , vertexColor_(assembler.gridGeometry().vertexMapper().size(), EntityColor::red)
112 {}
113
114 // returns number of green elements
115 std::size_t computeColors(const Assembler& assembler,
116 const std::vector<Scalar>& distanceFromLastLinearization,
117 Scalar threshold)
118 {
119 const auto& gridGeometry = assembler.gridGeometry();
120 const auto& gridView = gridGeometry.gridView();
121 const auto& elementMapper = gridGeometry.elementMapper();
122 const auto& vertexMapper = gridGeometry.vertexMapper();
123
124 // set all vertices to green
125 vertexColor_.assign(vertexColor_.size(), EntityColor::green);
126
127 // mark the red vertices
128 for (unsigned int i = 0; i < vertexColor_.size(); ++i)
129 {
130 using std::max;
131 if (distanceFromLastLinearization[i] > threshold)
132 // mark vertex as red if discrepancy is larger than
133 // the relative tolerance
134 vertexColor_[i] = EntityColor::red;
135 }
136
137 // Mark all red elements
138 for (const auto& element : elements(gridView))
139 {
140 // find out whether the current element features a red vertex
141 bool isRed = false;
142
143 int numVertices = element.subEntities(dim);
144
145 for (int i = 0; i < numVertices; ++i) {
146 int globalI = vertexMapper.subIndex(element, i, dim);
147
148 if (vertexColor_[globalI] == EntityColor::red) {
149 isRed = true;
150 break;
151 }
152 }
153
154 int eIdx = elementMapper.index(element);
155 // if a vertex is red, the element color is also red, otherwise green
156 if (isRed)
157 elementColor_[eIdx] = EntityColor::red;
158 else
159 elementColor_[eIdx] = EntityColor::green;
160 }
161
162 // mark orange vertices
163 for (const auto& element : elements(gridView))
164 {
165 int eIdx = elementMapper.index(element);
166
167 // only red elements tint vertices yellow
168 if (elementColor_[eIdx] == EntityColor::red)
169 {
170 int numVertices = element.subEntities(dim);
171
172 for (int i = 0; i < numVertices; ++i) {
173 int globalI = vertexMapper.subIndex(element, i, dim);
174
175 // red vertices don't become orange
176 if (vertexColor_[globalI] != EntityColor::red)
177 vertexColor_[globalI] = EntityColor::orange;
178 }
179 }
180 }
181
182 // at this point we communicate the yellow vertices to the
183 // neighboring processes because a neigbor process may not see
184 // the red vertex for yellow border vertices
186 minHandle(vertexMapper, vertexColor_);
187 gridView.communicate(minHandle,
188 Dune::InteriorBorder_InteriorBorder_Interface,
189 Dune::ForwardCommunication);
190
191 // mark yellow elements
192 for (const auto& element : elements(gridView))
193 {
194 int eIdx = elementMapper.index(element);
195
196 // only treat non-red elements
197 if (elementColor_[eIdx] != EntityColor::red)
198 {
199 // check whether the element features a orange vertex
200 bool isOrange = false;
201 int numVertices = element.subEntities(dim);
202
203 for (int i = 0; i < numVertices; ++i) {
204 int globalI = vertexMapper.subIndex(element, i, dim);
205
206 if (vertexColor_[globalI] == EntityColor::orange) {
207 isOrange = true;
208 break;
209 }
210 }
211
212 if (isOrange)
213 elementColor_[eIdx] = EntityColor::yellow;
214 }
215 }
216
217 // change orange vertices to yellow ones if it has at least
218 // one green element as a neighbor
219 for (const auto& element : elements(gridView))
220 {
221 int eIdx = elementMapper.index(element);
222
223 // only green elements are considered
224 if (elementColor_[eIdx] == EntityColor::green)
225 {
226 int numVertices = element.subEntities(dim);
227
228 for (int i = 0; i < numVertices; ++i) {
229 int globalI = vertexMapper.subIndex(element, i, dim);
230
231 // if a vertex is orange, recolor it to yellow
232 if (vertexColor_[globalI] == EntityColor::orange)
233 vertexColor_[globalI] = EntityColor::yellow;
234 }
235 }
236 }
237
238 // demote the border orange vertices
240 maxHandle(vertexMapper, vertexColor_);
241 gridView.communicate(maxHandle,
242 Dune::InteriorBorder_InteriorBorder_Interface,
243 Dune::ForwardCommunication);
244
245 // promote the remaining orange vertices to red
246 for (unsigned int i=0; i < vertexColor_.size(); ++i) {
247 // if a vertex is green or yellow don't do anything
248 if (vertexColor_[i] == EntityColor::green || vertexColor_[i] == EntityColor::yellow)
249 continue;
250
251 // set the vertex to red
252 vertexColor_[i] = EntityColor::red;
253 }
254
255 // count green elements
256 return std::count_if(elementColor_.begin(), elementColor_.end(),
257 [](EntityColor c){ return c == EntityColor::green; });
258 }
259
260 void resetJacobian(Assembler& assembler) const
261 {
262 auto& jacobian = assembler.jacobian();
263
264 // loop over all dofs
265 for (unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
266 {
267 // reset all entries corrosponding to a non-green vertex
268 if (vertexColor_[rowIdx] != EntityColor::green)
269 {
270 // set all matrix entries in the row to 0
271 auto colIt = jacobian[rowIdx].begin();
272 const auto& colEndIt = jacobian[rowIdx].end();
273 for (; colIt != colEndIt; ++colIt) {
274 *colIt = 0.0;
275 }
276 }
277 }
278 }
279
281 {
282 elementColor_.assign(elementColor_.size(), EntityColor::red);
283 vertexColor_.assign(vertexColor_.size(), EntityColor::red);
284 }
285
286 EntityColor elementColor(size_t idx) const
287 { return elementColor_[idx]; }
288
289 EntityColor vertexColor(size_t idx) const
290 { return vertexColor_[idx]; }
291
292 EntityColor dofColor(size_t idx) const
293 { return vertexColor_[idx]; }
294
295private:
297 std::vector<EntityColor> elementColor_;
298 std::vector<EntityColor> vertexColor_;
299};
300
305template<class Assembler>
307{
308 using Scalar = typename Assembler::Scalar;
309 using GridGeometry = typename Assembler::GridGeometry;
310 using JacobianMatrix = typename Assembler::JacobianMatrix;
311
312public:
313 PartialReassemblerEngine(const Assembler& assembler)
314 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
315 {}
316
317 // returns number of green elements
318 std::size_t computeColors(const Assembler& assembler,
319 const std::vector<Scalar>& distanceFromLastLinearization,
320 Scalar threshold)
321 {
322 const auto& gridGeometry = assembler.gridGeometry();
323 const auto& gridView = gridGeometry.gridView();
324 const auto& elementMapper = gridGeometry.elementMapper();
325
326 // mark the red elements
327 for (const auto& element : elements(gridView))
328 {
329 int eIdx = elementMapper.index(element);
330
331 if (distanceFromLastLinearization[eIdx] > threshold)
332 {
333 // mark element as red if discrepancy is larger than
334 // the relative tolerance
335 elementColor_[eIdx] = EntityColor::red;
336 }
337 else
338 {
339 elementColor_[eIdx] = EntityColor::green;
340 }
341 }
342
343 // mark the neighbors also red
344 const auto& connectivityMap = gridGeometry.connectivityMap();
345 for (unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
346 {
347 if (elementColor_[eIdx] == EntityColor::red)
348 continue; // element is red already!
349
350 if (distanceFromLastLinearization[eIdx] > threshold)
351 {
352 for (const auto& connectedDof : connectivityMap[eIdx])
353 elementColor_[connectedDof.globalJ] = EntityColor::red;
354 }
355 }
356
357 // count green elements
358 return std::count_if(elementColor_.begin(), elementColor_.end(),
359 [](EntityColor c){return c == EntityColor::green;});
360
361 }
362
363 void resetJacobian(Assembler& assembler) const
364 {
365 auto& jacobian = assembler.jacobian();
366 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
367
368 // loop over all dofs
369 for (unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
370 {
371 // reset all entries corresponding to a non-green element
372 if (elementColor_[colIdx] != EntityColor::green)
373 {
374 // set all matrix entries in the column to 0
375 jacobian[colIdx][colIdx] = 0;
376 for (const auto& dataJ : connectivityMap[colIdx])
377 jacobian[dataJ.globalJ][colIdx] = 0;
378 }
379 }
380 }
381
383 {
384 elementColor_.assign(elementColor_.size(), EntityColor::red);
385 }
386
387 EntityColor elementColor(size_t idx) const
388 { return elementColor_[idx]; }
389
390 EntityColor dofColor(size_t idx) const
391 { return elementColor_[idx]; }
392
393private:
395 std::vector<EntityColor> elementColor_;
396};
397
402template<class Assembler>
404: public PartialReassemblerEngine<Assembler, DiscretizationMethod::cctpfa>
405{
407public:
408 using ParentType::ParentType;
409};
410
413{
414 template<class Engine>
415 auto operator()(Engine&& e) -> decltype(e.vertexColor(0)) {}
416};
417
423template<class Assembler>
425{
426 using Scalar = typename Assembler::Scalar;
427 using GridGeometry = typename Assembler::GridGeometry;
428 using JacobianMatrix = typename Assembler::JacobianMatrix;
429 using VertexMapper = typename GridGeometry::VertexMapper;
430
431 static constexpr DiscretizationMethod discMethod = GridGeometry::discMethod;
433
434public:
435
440 PartialReassembler(const Assembler& assembler)
441 : engine_(assembler)
442 , greenElems_(0)
443 {
444 const auto& gridGeometry = assembler.gridGeometry();
445 totalElems_ = gridGeometry.elementMapper().size();
446 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
447 }
448
464 void computeColors(const Assembler& assembler,
465 const std::vector<Scalar>& distanceFromLastLinearization,
466 Scalar threshold)
467 {
468 greenElems_ = engine_.computeColors(assembler, distanceFromLastLinearization, threshold);
469 }
470
472 {
473 engine_.resetColors();
474 }
475
476 void resetJacobian(Assembler& assembler) const
477 {
478 engine_.resetJacobian(assembler);
479 }
480
484 template <class Communication>
485 void report(const Communication& comm, std::ostream& outStream)
486 {
487 if (comm.size() > 1)
488 greenElems_ = comm.sum(greenElems_);
489
490 const auto reassembledElems = totalElems_ - greenElems_;
491 const auto width = Fmt::formatted_size("{}", totalElems_);
492 outStream << Fmt::format(", reassembled {:{}} ({:3}%) elements",
493 reassembledElems, width, 100*reassembledElems/totalElems_);
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.
Formatting based on the fmt-library which implements std::format of C++20.
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:43
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:61
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:53
void computeColors(Args &&... args)
Definition: partialreassembler.hh:56
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:67
void report(Args &&... args)
Definition: partialreassembler.hh:50
DefaultPartialReassembler(Args &&... args)
Definition: partialreassembler.hh:46
void resetColors(Args &&... args)
Definition: partialreassembler.hh:59
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:64
the partial reassembler engine specialized for discretization methods
Definition: partialreassembler.hh:74
std::size_t computeColors(Args &&... args)
Definition: partialreassembler.hh:86
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:79
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:82
void resetColors(Args &&... args)
Definition: partialreassembler.hh:92
PartialReassemblerEngine(const Assembler &)
Definition: partialreassembler.hh:76
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:89
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:292
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:109
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:289
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:286
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:115
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:260
The partial reassembler engine specialized for the cellcentered TPFA method.
Definition: partialreassembler.hh:307
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:363
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:387
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:313
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:318
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:390
helper struct to determine whether the an engine class has vertex colors
Definition: partialreassembler.hh:413
auto operator()(Engine &&e) -> decltype(e.vertexColor(0))
Definition: partialreassembler.hh:415
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:425
PartialReassembler(const Assembler &assembler)
constructor
Definition: partialreassembler.hh:440
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:504
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:476
void resetColors()
Definition: partialreassembler.hh:471
void computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Determine the colors of entities for partial reassembly.
Definition: partialreassembler.hh:464
void report(const Communication &comm, std::ostream &outStream)
called by the assembler after successful assembly
Definition: partialreassembler.hh:485
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:78