24#ifndef DUMUX_PARTIAL_REASSEMBLER_HH
25#define DUMUX_PARTIAL_REASSEMBLER_HH
30#include <dune/grid/common/gridenums.hh>
31#include <dune/istl/multitypeblockvector.hh>
45 template<
typename... Args>
47 { DUNE_THROW(Dune::InvalidStateException,
"DefaultPartialReassembler should be never constructed!"); }
49 template<
typename... Args>
52 template<
typename... Args>
55 template<
typename... Args>
58 template<
typename... Args>
72template<
class Assembler, DiscretizationMethod discMethod>
77 { DUNE_THROW(Dune::NotImplemented,
"PartialReassembler for this discretization method!"); }
85 template<
typename... Args>
88 template<
typename... Args>
91 template<
typename... Args>
99template<
class Assembler>
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;
110 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
111 , vertexColor_(assembler.gridGeometry().vertexMapper().size(),
EntityColor::
red)
116 const std::vector<Scalar>& distanceFromLastLinearization,
119 const auto& gridGeometry = assembler.gridGeometry();
120 const auto& gridView = gridGeometry.gridView();
121 const auto& elementMapper = gridGeometry.elementMapper();
122 const auto& vertexMapper = gridGeometry.vertexMapper();
128 for (
unsigned int i = 0; i < vertexColor_.size(); ++i)
131 if (distanceFromLastLinearization[i] > threshold)
138 for (
const auto& element : elements(gridView))
143 int numVertices = element.subEntities(dim);
145 for (
int i = 0; i < numVertices; ++i) {
146 int globalI = vertexMapper.subIndex(element, i, dim);
154 int eIdx = elementMapper.index(element);
163 for (
const auto& element : elements(gridView))
165 int eIdx = elementMapper.index(element);
170 int numVertices = element.subEntities(dim);
172 for (
int i = 0; i < numVertices; ++i) {
173 int globalI = vertexMapper.subIndex(element, i, dim);
186 minHandle(vertexMapper, vertexColor_);
187 gridView.communicate(minHandle,
188 Dune::InteriorBorder_InteriorBorder_Interface,
189 Dune::ForwardCommunication);
192 for (
const auto& element : elements(gridView))
194 int eIdx = elementMapper.index(element);
200 bool isOrange =
false;
201 int numVertices = element.subEntities(dim);
203 for (
int i = 0; i < numVertices; ++i) {
204 int globalI = vertexMapper.subIndex(element, i, dim);
219 for (
const auto& element : elements(gridView))
221 int eIdx = elementMapper.index(element);
226 int numVertices = element.subEntities(dim);
228 for (
int i = 0; i < numVertices; ++i) {
229 int globalI = vertexMapper.subIndex(element, i, dim);
240 maxHandle(vertexMapper, vertexColor_);
241 gridView.communicate(maxHandle,
242 Dune::InteriorBorder_InteriorBorder_Interface,
243 Dune::ForwardCommunication);
246 for (
unsigned int i=0; i < vertexColor_.size(); ++i) {
256 return std::count_if(elementColor_.begin(), elementColor_.end(),
257 [](
EntityColor c){ return c == EntityColor::green; });
262 auto& jacobian = assembler.jacobian();
265 for (
unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
271 auto colIt = jacobian[rowIdx].begin();
272 const auto& colEndIt = jacobian[rowIdx].end();
273 for (; colIt != colEndIt; ++colIt) {
287 {
return elementColor_[idx]; }
290 {
return vertexColor_[idx]; }
293 {
return vertexColor_[idx]; }
297 std::vector<EntityColor> elementColor_;
298 std::vector<EntityColor> vertexColor_;
305template<
class Assembler>
308 using Scalar =
typename Assembler::Scalar;
309 using GridGeometry =
typename Assembler::GridGeometry;
310 using JacobianMatrix =
typename Assembler::JacobianMatrix;
314 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
319 const std::vector<Scalar>& distanceFromLastLinearization,
322 const auto& gridGeometry = assembler.gridGeometry();
323 const auto& gridView = gridGeometry.gridView();
324 const auto& elementMapper = gridGeometry.elementMapper();
327 for (
const auto& element : elements(gridView))
329 int eIdx = elementMapper.index(element);
331 if (distanceFromLastLinearization[eIdx] > threshold)
344 const auto& connectivityMap = gridGeometry.connectivityMap();
345 for (
unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
350 if (distanceFromLastLinearization[eIdx] > threshold)
352 for (
const auto& connectedDof : connectivityMap[eIdx])
358 return std::count_if(elementColor_.begin(), elementColor_.end(),
359 [](
EntityColor c){return c == EntityColor::green;});
365 auto& jacobian = assembler.jacobian();
366 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
369 for (
unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
375 jacobian[colIdx][colIdx] = 0;
376 for (
const auto& dataJ : connectivityMap[colIdx])
377 jacobian[dataJ.globalJ][colIdx] = 0;
388 {
return elementColor_[idx]; }
391 {
return elementColor_[idx]; }
395 std::vector<EntityColor> elementColor_;
402template<
class Assembler>
408 using ParentType::ParentType;
414 template<
class Engine>
423template<
class Assembler>
426 using Scalar =
typename Assembler::Scalar;
427 using GridGeometry =
typename Assembler::GridGeometry;
428 using JacobianMatrix =
typename Assembler::JacobianMatrix;
429 using VertexMapper =
typename GridGeometry::VertexMapper;
444 const auto& gridGeometry = assembler.gridGeometry();
445 totalElems_ = gridGeometry.elementMapper().size();
446 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
465 const std::vector<Scalar>& distanceFromLastLinearization,
468 greenElems_ = engine_.
computeColors(assembler, distanceFromLastLinearization, threshold);
484 template <
class Communication>
485 void report(
const Communication& comm, std::ostream& outStream)
488 greenElems_ = comm.sum(greenElems_);
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_);
503 typename std::enable_if_t<enable, int> = 0>
505 {
return engine_.vertexColor(idx); }
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
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
void resetColors()
Definition: partialreassembler.hh:280
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
void resetColors()
Definition: partialreassembler.hh:382
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