version 3.10-dev
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_PARTIAL_REASSEMBLER_HH
13#define DUMUX_PARTIAL_REASSEMBLER_HH
14
15#include <algorithm>
16#include <vector>
17
18#include <dune/grid/common/gridenums.hh>
19#include <dune/istl/multitypeblockvector.hh>
20
21#include <dumux/io/format.hh>
26
27#include "entitycolor.hh"
28
29namespace Dumux {
30
32{
33public:
34 template<typename... Args>
36 { DUNE_THROW(Dune::InvalidStateException, "DefaultPartialReassembler should be never constructed!"); }
37
38 template<typename... Args>
39 void report(Args&&... args) {}
40
41 template<typename... Args>
42 void resetJacobian(Args&&... args) const {}
43
44 template<typename... Args>
45 void computeColors(Args&&... args) {}
46
47 template<typename... Args>
48 void resetColors(Args&&... args) {}
49
50 EntityColor dofColor(size_t idx) const
51 { return EntityColor::red; }
52
53 EntityColor elementColor(size_t idx) const
54 { return EntityColor::red; }
55
56 EntityColor vertexColor(size_t idx) const
57 { return EntityColor::red; }
58};
59
61template<class Assembler, class DiscretizationMethod>
63{
64public:
65 PartialReassemblerEngine(const Assembler&)
66 { DUNE_THROW(Dune::NotImplemented, "PartialReassembler for this discretization method!"); }
67
68 EntityColor elementColor(size_t idx) const
69 { return EntityColor::red; }
70
71 EntityColor dofColor(size_t idx) const
72 { return EntityColor::red; }
73
74 template<typename... Args>
75 std::size_t computeColors(Args&&... args) { return 0; }
76
77 template<typename... Args>
78 void resetJacobian(Args&&... args) const {}
79
80 template<typename... Args>
81 void resetColors(Args&&... args) {}
82};
83
88template<class Assembler>
89class PartialReassemblerEngine<Assembler, DiscretizationMethods::Box>
90{
91 using Scalar = typename Assembler::Scalar;
92 using GridGeometry = typename Assembler::GridGeometry;
93 using JacobianMatrix = typename Assembler::JacobianMatrix;
94 using VertexMapper = typename GridGeometry::VertexMapper;
95 static constexpr int dim = GridGeometry::GridView::dimension;
96
97public:
98 PartialReassemblerEngine(const Assembler& assembler)
99 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
100 , vertexColor_(assembler.gridGeometry().vertexMapper().size(), EntityColor::red)
101 {}
102
103 // returns number of green elements
104 std::size_t computeColors(const Assembler& assembler,
105 const std::vector<Scalar>& distanceFromLastLinearization,
106 Scalar threshold)
107 {
108 const auto& gridGeometry = assembler.gridGeometry();
109 const auto& gridView = gridGeometry.gridView();
110 const auto& elementMapper = gridGeometry.elementMapper();
111 const auto& vertexMapper = gridGeometry.vertexMapper();
112
113 // set all vertices to green
114 vertexColor_.assign(vertexColor_.size(), EntityColor::green);
115
116 // mark the red vertices
117 for (unsigned int i = 0; i < vertexColor_.size(); ++i)
118 {
119 using std::max;
120 if (distanceFromLastLinearization[i] > threshold)
121 // mark vertex as red if discrepancy is larger than
122 // the relative tolerance
123 vertexColor_[i] = EntityColor::red;
124 }
125
126 // Mark all red elements
127 for (const auto& element : elements(gridView))
128 {
129 // find out whether the current element features a red vertex
130 bool isRed = false;
131
132 int numVertices = element.subEntities(dim);
133
134 for (int i = 0; i < numVertices; ++i) {
135 int globalI = vertexMapper.subIndex(element, i, dim);
136
137 if (vertexColor_[globalI] == EntityColor::red) {
138 isRed = true;
139 break;
140 }
141 }
142
143 int eIdx = elementMapper.index(element);
144 // if a vertex is red, the element color is also red, otherwise green
145 if (isRed)
146 elementColor_[eIdx] = EntityColor::red;
147 else
148 elementColor_[eIdx] = EntityColor::green;
149 }
150
151 // mark orange vertices
152 for (const auto& element : elements(gridView))
153 {
154 int eIdx = elementMapper.index(element);
155
156 // only red elements tint vertices yellow
157 if (elementColor_[eIdx] == EntityColor::red)
158 {
159 int numVertices = element.subEntities(dim);
160
161 for (int i = 0; i < numVertices; ++i) {
162 int globalI = vertexMapper.subIndex(element, i, dim);
163
164 // red vertices don't become orange
165 if (vertexColor_[globalI] != EntityColor::red)
166 vertexColor_[globalI] = EntityColor::orange;
167 }
168 }
169 }
170
171 // at this point we communicate the yellow vertices to the
172 // neighboring processes because a neighbor process may not see
173 // the red vertex for yellow border vertices
175 minHandle(vertexMapper, vertexColor_);
176 if constexpr (Detail::canCommunicate<typename GridGeometry::GridView::Traits::Grid, dim>)
177 gridView.communicate(minHandle,
178 Dune::InteriorBorder_InteriorBorder_Interface,
179 Dune::ForwardCommunication);
180 else
181 DUNE_THROW(Dune::InvalidStateException, "Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim << "-entities.");
182
183 // mark yellow elements
184 for (const auto& element : elements(gridView))
185 {
186 int eIdx = elementMapper.index(element);
187
188 // only treat non-red elements
189 if (elementColor_[eIdx] != EntityColor::red)
190 {
191 // check whether the element features a orange vertex
192 bool isOrange = false;
193 int numVertices = element.subEntities(dim);
194
195 for (int i = 0; i < numVertices; ++i) {
196 int globalI = vertexMapper.subIndex(element, i, dim);
197
198 if (vertexColor_[globalI] == EntityColor::orange) {
199 isOrange = true;
200 break;
201 }
202 }
203
204 if (isOrange)
205 elementColor_[eIdx] = EntityColor::yellow;
206 }
207 }
208
209 // change orange vertices to yellow ones if it has at least
210 // one green element as a neighbor
211 for (const auto& element : elements(gridView))
212 {
213 int eIdx = elementMapper.index(element);
214
215 // only green elements are considered
216 if (elementColor_[eIdx] == EntityColor::green)
217 {
218 int numVertices = element.subEntities(dim);
219
220 for (int i = 0; i < numVertices; ++i) {
221 int globalI = vertexMapper.subIndex(element, i, dim);
222
223 // if a vertex is orange, recolor it to yellow
224 if (vertexColor_[globalI] == EntityColor::orange)
225 vertexColor_[globalI] = EntityColor::yellow;
226 }
227 }
228 }
229
230 // demote the border orange vertices
232 maxHandle(vertexMapper, vertexColor_);
233 if constexpr (Detail::canCommunicate<typename GridGeometry::GridView::Traits::Grid, dim>)
234 gridView.communicate(maxHandle,
235 Dune::InteriorBorder_InteriorBorder_Interface,
236 Dune::ForwardCommunication);
237 else
238 DUNE_THROW(Dune::InvalidStateException, "Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim << "-entities.");
239
240 // promote the remaining orange vertices to red
241 for (unsigned int i=0; i < vertexColor_.size(); ++i) {
242 // if a vertex is green or yellow don't do anything
243 if (vertexColor_[i] == EntityColor::green || vertexColor_[i] == EntityColor::yellow)
244 continue;
245
246 // set the vertex to red
247 vertexColor_[i] = EntityColor::red;
248 }
249
250 // count green elements
251 return std::count_if(elementColor_.begin(), elementColor_.end(),
252 [](EntityColor c){ return c == EntityColor::green; });
253 }
254
255 void resetJacobian(Assembler& assembler) const
256 {
257 auto& jacobian = assembler.jacobian();
258
259 // loop over all dofs
260 for (unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
261 {
262 // reset all entries corresponding to a non-green vertex
263 if (vertexColor_[rowIdx] != EntityColor::green)
264 {
265 // set all matrix entries in the row to 0
266 auto colIt = jacobian[rowIdx].begin();
267 const auto& colEndIt = jacobian[rowIdx].end();
268 for (; colIt != colEndIt; ++colIt) {
269 *colIt = 0.0;
270 }
271 }
272 }
273 }
274
276 {
277 elementColor_.assign(elementColor_.size(), EntityColor::red);
278 vertexColor_.assign(vertexColor_.size(), EntityColor::red);
279 }
280
281 EntityColor elementColor(size_t idx) const
282 { return elementColor_[idx]; }
283
284 EntityColor vertexColor(size_t idx) const
285 { return vertexColor_[idx]; }
286
287 EntityColor dofColor(size_t idx) const
288 { return vertexColor_[idx]; }
289
290private:
292 std::vector<EntityColor> elementColor_;
293 std::vector<EntityColor> vertexColor_;
294};
295
300template<class Assembler>
301class PartialReassemblerEngine<Assembler, DiscretizationMethods::CCTpfa>
302{
303 using Scalar = typename Assembler::Scalar;
304 using GridGeometry = typename Assembler::GridGeometry;
305 using JacobianMatrix = typename Assembler::JacobianMatrix;
306
307public:
308 PartialReassemblerEngine(const Assembler& assembler)
309 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
310 {}
311
312 // returns number of green elements
313 std::size_t computeColors(const Assembler& assembler,
314 const std::vector<Scalar>& distanceFromLastLinearization,
315 Scalar threshold)
316 {
317 const auto& gridGeometry = assembler.gridGeometry();
318 const auto& gridView = gridGeometry.gridView();
319 const auto& elementMapper = gridGeometry.elementMapper();
320
321 // mark the red elements
322 for (const auto& element : elements(gridView))
323 {
324 int eIdx = elementMapper.index(element);
325
326 if (distanceFromLastLinearization[eIdx] > threshold)
327 {
328 // mark element as red if discrepancy is larger than
329 // the relative tolerance
330 elementColor_[eIdx] = EntityColor::red;
331 }
332 else
333 {
334 elementColor_[eIdx] = EntityColor::green;
335 }
336 }
337
338 // mark the neighbors also red
339 const auto& connectivityMap = gridGeometry.connectivityMap();
340 for (unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
341 {
342 if (elementColor_[eIdx] == EntityColor::red)
343 continue; // element is red already!
344
345 if (distanceFromLastLinearization[eIdx] > threshold)
346 {
347 for (const auto& connectedDof : connectivityMap[eIdx])
348 elementColor_[connectedDof.globalJ] = EntityColor::red;
349 }
350 }
351
352 // count green elements
353 return std::count_if(elementColor_.begin(), elementColor_.end(),
354 [](EntityColor c){return c == EntityColor::green;});
355
356 }
357
358 void resetJacobian(Assembler& assembler) const
359 {
360 auto& jacobian = assembler.jacobian();
361 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
362
363 // loop over all dofs
364 for (unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
365 {
366 // reset all entries corresponding to a non-green element
367 if (elementColor_[colIdx] != EntityColor::green)
368 {
369 // set all matrix entries in the column to 0
370 jacobian[colIdx][colIdx] = 0;
371 for (const auto& dataJ : connectivityMap[colIdx])
372 jacobian[dataJ.globalJ][colIdx] = 0;
373 }
374 }
375 }
376
378 {
379 elementColor_.assign(elementColor_.size(), EntityColor::red);
380 }
381
382 EntityColor elementColor(size_t idx) const
383 { return elementColor_[idx]; }
384
385 EntityColor dofColor(size_t idx) const
386 { return elementColor_[idx]; }
387
388private:
390 std::vector<EntityColor> elementColor_;
391};
392
397template<class Assembler>
398class PartialReassemblerEngine<Assembler, DiscretizationMethods::CCMpfa>
399: public PartialReassemblerEngine<Assembler, DiscretizationMethods::CCTpfa>
400{
402public:
403 using ParentType::ParentType;
404};
405
408{
409 template<class Engine>
410 auto operator()(Engine&& e) -> decltype(e.vertexColor(0)) {}
411};
412
418template<class Assembler>
420{
421 using Scalar = typename Assembler::Scalar;
422 using GridGeometry = typename Assembler::GridGeometry;
423 using JacobianMatrix = typename Assembler::JacobianMatrix;
424 using VertexMapper = typename GridGeometry::VertexMapper;
425
426 using DiscretizationMethod = typename GridGeometry::DiscretizationMethod;
428
429public:
430
435 PartialReassembler(const Assembler& assembler)
436 : engine_(assembler)
437 , greenElems_(0)
438 {
439 const auto& gridGeometry = assembler.gridGeometry();
440 totalElems_ = gridGeometry.elementMapper().size();
441 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
442 }
443
459 void computeColors(const Assembler& assembler,
460 const std::vector<Scalar>& distanceFromLastLinearization,
461 Scalar threshold)
462 {
463 greenElems_ = engine_.computeColors(assembler, distanceFromLastLinearization, threshold);
464 }
465
467 {
468 engine_.resetColors();
469 }
470
471 void resetJacobian(Assembler& assembler) const
472 {
473 engine_.resetJacobian(assembler);
474 }
475
479 template <class Communication>
480 void report(const Communication& comm, std::ostream& outStream)
481 {
482 if (comm.size() > 1)
483 greenElems_ = comm.sum(greenElems_);
484
485 const auto reassembledElems = totalElems_ - greenElems_;
486 const auto width = Fmt::formatted_size("{}", totalElems_);
487 outStream << Fmt::format(", reassembled {:{}} ({:3}%) elements",
488 reassembledElems, width, 100*reassembledElems/totalElems_);
489 }
490
491 EntityColor elementColor(size_t idx) const
492 { return engine_.elementColor(idx); }
493
494 EntityColor dofColor(size_t idx) const
495 { return engine_.dofColor(idx); }
496
497 template<bool enable = decltype(isValid(hasVertexColor()).template check<Engine>())::value,
498 typename std::enable_if_t<enable, int> = 0>
499 EntityColor vertexColor(size_t idx) const
500 { return engine_.vertexColor(idx); }
501
502private:
503 Engine engine_;
504 size_t totalElems_;
505 size_t greenElems_;
506};
507
508} // namespace Dumux
509
510#endif
Definition: partialreassembler.hh:32
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:50
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:42
void computeColors(Args &&... args)
Definition: partialreassembler.hh:45
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:56
void report(Args &&... args)
Definition: partialreassembler.hh:39
DefaultPartialReassembler(Args &&... args)
Definition: partialreassembler.hh:35
void resetColors(Args &&... args)
Definition: partialreassembler.hh:48
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:53
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:287
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:104
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:98
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:255
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:284
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:281
The partial reassembler engine specialized for the cellcentered TPFA method.
Definition: partialreassembler.hh:302
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:358
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:313
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:385
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:308
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:382
the partial reassembler engine specialized for discretization methods
Definition: partialreassembler.hh:63
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:78
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:68
PartialReassemblerEngine(const Assembler &)
Definition: partialreassembler.hh:65
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:71
void resetColors(Args &&... args)
Definition: partialreassembler.hh:81
std::size_t computeColors(Args &&... args)
Definition: partialreassembler.hh:75
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
PartialReassembler(const Assembler &assembler)
constructor
Definition: partialreassembler.hh:435
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:499
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:471
void resetColors()
Definition: partialreassembler.hh:466
void computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Determine the colors of entities for partial reassembly.
Definition: partialreassembler.hh:459
void report(const Communication &comm, std::ostream &outStream)
called by the assembler after successful assembly
Definition: partialreassembler.hh:480
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:494
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:491
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:66
An enum class to define the colors of elements and vertices required for partial Jacobian reassembly.
Formatting based on the fmt-library which implements std::format of C++20.
dune-grid capabilities compatibility layer
EntityColor
The colors of elements and vertices required for partial Jacobian reassembly.
Definition: entitycolor.hh:23
@ green
does not need to be reassembled
@ red
distance from last linearization is above the tolerance
@ yellow
neighboring entity is red
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
Definition: adapt.hh:17
helper struct to determine whether the an engine class has vertex colors
Definition: partialreassembler.hh:408
auto operator()(Engine &&e) -> decltype(e.vertexColor(0))
Definition: partialreassembler.hh:410
Contains a class to exchange entries of a vector.