3.6-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
33#include <dumux/io/format.hh>
38
39#include "entitycolor.hh"
40
41namespace Dumux {
42
44{
45public:
46 template<typename... Args>
48 { DUNE_THROW(Dune::InvalidStateException, "DefaultPartialReassembler should be never constructed!"); }
49
50 template<typename... Args>
51 void report(Args&&... args) {}
52
53 template<typename... Args>
54 void resetJacobian(Args&&... args) const {}
55
56 template<typename... Args>
57 void computeColors(Args&&... args) {}
58
59 template<typename... Args>
60 void resetColors(Args&&... args) {}
61
62 EntityColor dofColor(size_t idx) const
63 { return EntityColor::red; }
64
65 EntityColor elementColor(size_t idx) const
66 { return EntityColor::red; }
67
68 EntityColor vertexColor(size_t idx) const
69 { return EntityColor::red; }
70};
71
73template<class Assembler, class DiscretizationMethod>
75{
76public:
77 PartialReassemblerEngine(const Assembler&)
78 { DUNE_THROW(Dune::NotImplemented, "PartialReassembler for this discretization method!"); }
79
80 EntityColor elementColor(size_t idx) const
81 { return EntityColor::red; }
82
83 EntityColor dofColor(size_t idx) const
84 { return EntityColor::red; }
85
86 template<typename... Args>
87 std::size_t computeColors(Args&&... args) { return 0; }
88
89 template<typename... Args>
90 void resetJacobian(Args&&... args) const {}
91
92 template<typename... Args>
93 void resetColors(Args&&... args) {}
94};
95
100template<class Assembler>
101class PartialReassemblerEngine<Assembler, DiscretizationMethods::Box>
102{
103 using Scalar = typename Assembler::Scalar;
104 using GridGeometry = typename Assembler::GridGeometry;
105 using JacobianMatrix = typename Assembler::JacobianMatrix;
106 using VertexMapper = typename GridGeometry::VertexMapper;
107 static constexpr int dim = GridGeometry::GridView::dimension;
108
109public:
110 PartialReassemblerEngine(const Assembler& assembler)
111 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
112 , vertexColor_(assembler.gridGeometry().vertexMapper().size(), EntityColor::red)
113 {}
114
115 // returns number of green elements
116 std::size_t computeColors(const Assembler& assembler,
117 const std::vector<Scalar>& distanceFromLastLinearization,
118 Scalar threshold)
119 {
120 const auto& gridGeometry = assembler.gridGeometry();
121 const auto& gridView = gridGeometry.gridView();
122 const auto& elementMapper = gridGeometry.elementMapper();
123 const auto& vertexMapper = gridGeometry.vertexMapper();
124
125 // set all vertices to green
126 vertexColor_.assign(vertexColor_.size(), EntityColor::green);
127
128 // mark the red vertices
129 for (unsigned int i = 0; i < vertexColor_.size(); ++i)
130 {
131 using std::max;
132 if (distanceFromLastLinearization[i] > threshold)
133 // mark vertex as red if discrepancy is larger than
134 // the relative tolerance
135 vertexColor_[i] = EntityColor::red;
136 }
137
138 // Mark all red elements
139 for (const auto& element : elements(gridView))
140 {
141 // find out whether the current element features a red vertex
142 bool isRed = false;
143
144 int numVertices = element.subEntities(dim);
145
146 for (int i = 0; i < numVertices; ++i) {
147 int globalI = vertexMapper.subIndex(element, i, dim);
148
149 if (vertexColor_[globalI] == EntityColor::red) {
150 isRed = true;
151 break;
152 }
153 }
154
155 int eIdx = elementMapper.index(element);
156 // if a vertex is red, the element color is also red, otherwise green
157 if (isRed)
158 elementColor_[eIdx] = EntityColor::red;
159 else
160 elementColor_[eIdx] = EntityColor::green;
161 }
162
163 // mark orange vertices
164 for (const auto& element : elements(gridView))
165 {
166 int eIdx = elementMapper.index(element);
167
168 // only red elements tint vertices yellow
169 if (elementColor_[eIdx] == EntityColor::red)
170 {
171 int numVertices = element.subEntities(dim);
172
173 for (int i = 0; i < numVertices; ++i) {
174 int globalI = vertexMapper.subIndex(element, i, dim);
175
176 // red vertices don't become orange
177 if (vertexColor_[globalI] != EntityColor::red)
178 vertexColor_[globalI] = EntityColor::orange;
179 }
180 }
181 }
182
183 // at this point we communicate the yellow vertices to the
184 // neighboring processes because a neighbor process may not see
185 // the red vertex for yellow border vertices
187 minHandle(vertexMapper, vertexColor_);
188 if constexpr (Detail::canCommunicate<typename GridGeometry::GridView::Traits::Grid, dim>)
189 gridView.communicate(minHandle,
190 Dune::InteriorBorder_InteriorBorder_Interface,
191 Dune::ForwardCommunication);
192 else
193 DUNE_THROW(Dune::InvalidStateException, "Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim << "-entities.");
194
195 // mark yellow elements
196 for (const auto& element : elements(gridView))
197 {
198 int eIdx = elementMapper.index(element);
199
200 // only treat non-red elements
201 if (elementColor_[eIdx] != EntityColor::red)
202 {
203 // check whether the element features a orange vertex
204 bool isOrange = false;
205 int numVertices = element.subEntities(dim);
206
207 for (int i = 0; i < numVertices; ++i) {
208 int globalI = vertexMapper.subIndex(element, i, dim);
209
210 if (vertexColor_[globalI] == EntityColor::orange) {
211 isOrange = true;
212 break;
213 }
214 }
215
216 if (isOrange)
217 elementColor_[eIdx] = EntityColor::yellow;
218 }
219 }
220
221 // change orange vertices to yellow ones if it has at least
222 // one green element as a neighbor
223 for (const auto& element : elements(gridView))
224 {
225 int eIdx = elementMapper.index(element);
226
227 // only green elements are considered
228 if (elementColor_[eIdx] == EntityColor::green)
229 {
230 int numVertices = element.subEntities(dim);
231
232 for (int i = 0; i < numVertices; ++i) {
233 int globalI = vertexMapper.subIndex(element, i, dim);
234
235 // if a vertex is orange, recolor it to yellow
236 if (vertexColor_[globalI] == EntityColor::orange)
237 vertexColor_[globalI] = EntityColor::yellow;
238 }
239 }
240 }
241
242 // demote the border orange vertices
244 maxHandle(vertexMapper, vertexColor_);
245 if constexpr (Detail::canCommunicate<typename GridGeometry::GridView::Traits::Grid, dim>)
246 gridView.communicate(maxHandle,
247 Dune::InteriorBorder_InteriorBorder_Interface,
248 Dune::ForwardCommunication);
249 else
250 DUNE_THROW(Dune::InvalidStateException, "Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim << "-entities.");
251
252 // promote the remaining orange vertices to red
253 for (unsigned int i=0; i < vertexColor_.size(); ++i) {
254 // if a vertex is green or yellow don't do anything
255 if (vertexColor_[i] == EntityColor::green || vertexColor_[i] == EntityColor::yellow)
256 continue;
257
258 // set the vertex to red
259 vertexColor_[i] = EntityColor::red;
260 }
261
262 // count green elements
263 return std::count_if(elementColor_.begin(), elementColor_.end(),
264 [](EntityColor c){ return c == EntityColor::green; });
265 }
266
267 void resetJacobian(Assembler& assembler) const
268 {
269 auto& jacobian = assembler.jacobian();
270
271 // loop over all dofs
272 for (unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
273 {
274 // reset all entries corresponding to a non-green vertex
275 if (vertexColor_[rowIdx] != EntityColor::green)
276 {
277 // set all matrix entries in the row to 0
278 auto colIt = jacobian[rowIdx].begin();
279 const auto& colEndIt = jacobian[rowIdx].end();
280 for (; colIt != colEndIt; ++colIt) {
281 *colIt = 0.0;
282 }
283 }
284 }
285 }
286
288 {
289 elementColor_.assign(elementColor_.size(), EntityColor::red);
290 vertexColor_.assign(vertexColor_.size(), EntityColor::red);
291 }
292
293 EntityColor elementColor(size_t idx) const
294 { return elementColor_[idx]; }
295
296 EntityColor vertexColor(size_t idx) const
297 { return vertexColor_[idx]; }
298
299 EntityColor dofColor(size_t idx) const
300 { return vertexColor_[idx]; }
301
302private:
304 std::vector<EntityColor> elementColor_;
305 std::vector<EntityColor> vertexColor_;
306};
307
312template<class Assembler>
313class PartialReassemblerEngine<Assembler, DiscretizationMethods::CCTpfa>
314{
315 using Scalar = typename Assembler::Scalar;
316 using GridGeometry = typename Assembler::GridGeometry;
317 using JacobianMatrix = typename Assembler::JacobianMatrix;
318
319public:
320 PartialReassemblerEngine(const Assembler& assembler)
321 : elementColor_(assembler.gridGeometry().elementMapper().size(), EntityColor::red)
322 {}
323
324 // returns number of green elements
325 std::size_t computeColors(const Assembler& assembler,
326 const std::vector<Scalar>& distanceFromLastLinearization,
327 Scalar threshold)
328 {
329 const auto& gridGeometry = assembler.gridGeometry();
330 const auto& gridView = gridGeometry.gridView();
331 const auto& elementMapper = gridGeometry.elementMapper();
332
333 // mark the red elements
334 for (const auto& element : elements(gridView))
335 {
336 int eIdx = elementMapper.index(element);
337
338 if (distanceFromLastLinearization[eIdx] > threshold)
339 {
340 // mark element as red if discrepancy is larger than
341 // the relative tolerance
342 elementColor_[eIdx] = EntityColor::red;
343 }
344 else
345 {
346 elementColor_[eIdx] = EntityColor::green;
347 }
348 }
349
350 // mark the neighbors also red
351 const auto& connectivityMap = gridGeometry.connectivityMap();
352 for (unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
353 {
354 if (elementColor_[eIdx] == EntityColor::red)
355 continue; // element is red already!
356
357 if (distanceFromLastLinearization[eIdx] > threshold)
358 {
359 for (const auto& connectedDof : connectivityMap[eIdx])
360 elementColor_[connectedDof.globalJ] = EntityColor::red;
361 }
362 }
363
364 // count green elements
365 return std::count_if(elementColor_.begin(), elementColor_.end(),
366 [](EntityColor c){return c == EntityColor::green;});
367
368 }
369
370 void resetJacobian(Assembler& assembler) const
371 {
372 auto& jacobian = assembler.jacobian();
373 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
374
375 // loop over all dofs
376 for (unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
377 {
378 // reset all entries corresponding to a non-green element
379 if (elementColor_[colIdx] != EntityColor::green)
380 {
381 // set all matrix entries in the column to 0
382 jacobian[colIdx][colIdx] = 0;
383 for (const auto& dataJ : connectivityMap[colIdx])
384 jacobian[dataJ.globalJ][colIdx] = 0;
385 }
386 }
387 }
388
390 {
391 elementColor_.assign(elementColor_.size(), EntityColor::red);
392 }
393
394 EntityColor elementColor(size_t idx) const
395 { return elementColor_[idx]; }
396
397 EntityColor dofColor(size_t idx) const
398 { return elementColor_[idx]; }
399
400private:
402 std::vector<EntityColor> elementColor_;
403};
404
409template<class Assembler>
410class PartialReassemblerEngine<Assembler, DiscretizationMethods::CCMpfa>
411: public PartialReassemblerEngine<Assembler, DiscretizationMethods::CCTpfa>
412{
414public:
415 using ParentType::ParentType;
416};
417
420{
421 template<class Engine>
422 auto operator()(Engine&& e) -> decltype(e.vertexColor(0)) {}
423};
424
430template<class Assembler>
432{
433 using Scalar = typename Assembler::Scalar;
434 using GridGeometry = typename Assembler::GridGeometry;
435 using JacobianMatrix = typename Assembler::JacobianMatrix;
436 using VertexMapper = typename GridGeometry::VertexMapper;
437
438 using DiscretizationMethod = typename GridGeometry::DiscretizationMethod;
440
441public:
442
447 PartialReassembler(const Assembler& assembler)
448 : engine_(assembler)
449 , greenElems_(0)
450 {
451 const auto& gridGeometry = assembler.gridGeometry();
452 totalElems_ = gridGeometry.elementMapper().size();
453 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
454 }
455
471 void computeColors(const Assembler& assembler,
472 const std::vector<Scalar>& distanceFromLastLinearization,
473 Scalar threshold)
474 {
475 greenElems_ = engine_.computeColors(assembler, distanceFromLastLinearization, threshold);
476 }
477
479 {
480 engine_.resetColors();
481 }
482
483 void resetJacobian(Assembler& assembler) const
484 {
485 engine_.resetJacobian(assembler);
486 }
487
491 template <class Communication>
492 void report(const Communication& comm, std::ostream& outStream)
493 {
494 if (comm.size() > 1)
495 greenElems_ = comm.sum(greenElems_);
496
497 const auto reassembledElems = totalElems_ - greenElems_;
498 const auto width = Fmt::formatted_size("{}", totalElems_);
499 outStream << Fmt::format(", reassembled {:{}} ({:3}%) elements",
500 reassembledElems, width, 100*reassembledElems/totalElems_);
501 }
502
503 EntityColor elementColor(size_t idx) const
504 { return engine_.elementColor(idx); }
505
506 EntityColor dofColor(size_t idx) const
507 { return engine_.dofColor(idx); }
508
509 template<bool enable = decltype(isValid(hasVertexColor()).template check<Engine>())::value,
510 typename std::enable_if_t<enable, int> = 0>
511 EntityColor vertexColor(size_t idx) const
512 { return engine_.vertexColor(idx); }
513
514private:
515 Engine engine_;
516 size_t totalElems_;
517 size_t greenElems_;
518};
519
520} // namespace Dumux
521
522#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.
dune-grid capabilities compatibility layer
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.
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
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
Definition: partialreassembler.hh:44
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:62
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:54
void computeColors(Args &&... args)
Definition: partialreassembler.hh:57
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:68
void report(Args &&... args)
Definition: partialreassembler.hh:51
DefaultPartialReassembler(Args &&... args)
Definition: partialreassembler.hh:47
void resetColors(Args &&... args)
Definition: partialreassembler.hh:60
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:65
the partial reassembler engine specialized for discretization methods
Definition: partialreassembler.hh:75
void resetJacobian(Args &&... args) const
Definition: partialreassembler.hh:90
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:80
PartialReassemblerEngine(const Assembler &)
Definition: partialreassembler.hh:77
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:83
void resetColors(Args &&... args)
Definition: partialreassembler.hh:93
std::size_t computeColors(Args &&... args)
Definition: partialreassembler.hh:87
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:299
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:116
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:110
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:267
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:296
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:293
The partial reassembler engine specialized for the cellcentered TPFA method.
Definition: partialreassembler.hh:314
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:370
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition: partialreassembler.hh:325
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:397
PartialReassemblerEngine(const Assembler &assembler)
Definition: partialreassembler.hh:320
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:394
helper struct to determine whether the an engine class has vertex colors
Definition: partialreassembler.hh:420
auto operator()(Engine &&e) -> decltype(e.vertexColor(0))
Definition: partialreassembler.hh:422
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:432
PartialReassembler(const Assembler &assembler)
constructor
Definition: partialreassembler.hh:447
EntityColor vertexColor(size_t idx) const
Definition: partialreassembler.hh:511
void resetJacobian(Assembler &assembler) const
Definition: partialreassembler.hh:483
void resetColors()
Definition: partialreassembler.hh:478
void computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Determine the colors of entities for partial reassembly.
Definition: partialreassembler.hh:471
void report(const Communication &comm, std::ostream &outStream)
called by the assembler after successful assembly
Definition: partialreassembler.hh:492
EntityColor dofColor(size_t idx) const
Definition: partialreassembler.hh:506
EntityColor elementColor(size_t idx) const
Definition: partialreassembler.hh:503
A data handle class to exchange entries of a vector.
Definition: vectorcommdatahandle.hh:78