12#ifndef DUMUX_PARTIAL_REASSEMBLER_HH
13#define DUMUX_PARTIAL_REASSEMBLER_HH
18#include <dune/grid/common/gridenums.hh>
19#include <dune/istl/multitypeblockvector.hh>
34 template<
typename... Args>
36 { DUNE_THROW(Dune::InvalidStateException,
"DefaultPartialReassembler should be never constructed!"); }
38 template<
typename... Args>
41 template<
typename... Args>
44 template<
typename... Args>
47 template<
typename... Args>
61template<
class Assembler,
class DiscretizationMethod>
66 { DUNE_THROW(Dune::NotImplemented,
"PartialReassembler for this discretization method!"); }
74 template<
typename... Args>
77 template<
typename... Args>
80 template<
typename... Args>
88template<
class Assembler>
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;
99 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
100 , vertexColor_(assembler.gridGeometry().vertexMapper().size(),
EntityColor::
red)
105 const std::vector<Scalar>& distanceFromLastLinearization,
108 const auto& gridGeometry = assembler.gridGeometry();
109 const auto& gridView = gridGeometry.gridView();
110 const auto& elementMapper = gridGeometry.elementMapper();
111 const auto& vertexMapper = gridGeometry.vertexMapper();
117 for (
unsigned int i = 0; i < vertexColor_.size(); ++i)
120 if (distanceFromLastLinearization[i] > threshold)
127 for (
const auto& element : elements(gridView))
132 int numVertices = element.subEntities(dim);
134 for (
int i = 0; i < numVertices; ++i) {
135 int globalI = vertexMapper.subIndex(element, i, dim);
143 int eIdx = elementMapper.index(element);
152 for (
const auto& element : elements(gridView))
154 int eIdx = elementMapper.index(element);
159 int numVertices = element.subEntities(dim);
161 for (
int i = 0; i < numVertices; ++i) {
162 int globalI = vertexMapper.subIndex(element, i, dim);
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);
181 DUNE_THROW(Dune::InvalidStateException,
"Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim <<
"-entities.");
184 for (
const auto& element : elements(gridView))
186 int eIdx = elementMapper.index(element);
192 bool isOrange =
false;
193 int numVertices = element.subEntities(dim);
195 for (
int i = 0; i < numVertices; ++i) {
196 int globalI = vertexMapper.subIndex(element, i, dim);
211 for (
const auto& element : elements(gridView))
213 int eIdx = elementMapper.index(element);
218 int numVertices = element.subEntities(dim);
220 for (
int i = 0; i < numVertices; ++i) {
221 int globalI = vertexMapper.subIndex(element, i, dim);
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);
238 DUNE_THROW(Dune::InvalidStateException,
"Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim <<
"-entities.");
241 for (
unsigned int i=0; i < vertexColor_.size(); ++i) {
251 return std::count_if(elementColor_.begin(), elementColor_.end(),
252 [](
EntityColor c){ return c == EntityColor::green; });
257 auto& jacobian = assembler.jacobian();
260 for (
unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
266 auto colIt = jacobian[rowIdx].begin();
267 const auto& colEndIt = jacobian[rowIdx].end();
268 for (; colIt != colEndIt; ++colIt) {
282 {
return elementColor_[idx]; }
285 {
return vertexColor_[idx]; }
288 {
return vertexColor_[idx]; }
292 std::vector<EntityColor> elementColor_;
293 std::vector<EntityColor> vertexColor_;
300template<
class Assembler>
303 using Scalar =
typename Assembler::Scalar;
304 using GridGeometry =
typename Assembler::GridGeometry;
305 using JacobianMatrix =
typename Assembler::JacobianMatrix;
309 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
314 const std::vector<Scalar>& distanceFromLastLinearization,
317 const auto& gridGeometry = assembler.gridGeometry();
318 const auto& gridView = gridGeometry.gridView();
319 const auto& elementMapper = gridGeometry.elementMapper();
322 for (
const auto& element : elements(gridView))
324 int eIdx = elementMapper.index(element);
326 if (distanceFromLastLinearization[eIdx] > threshold)
339 const auto& connectivityMap = gridGeometry.connectivityMap();
340 for (
unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
345 if (distanceFromLastLinearization[eIdx] > threshold)
347 for (
const auto& connectedDof : connectivityMap[eIdx])
353 return std::count_if(elementColor_.begin(), elementColor_.end(),
354 [](
EntityColor c){return c == EntityColor::green;});
360 auto& jacobian = assembler.jacobian();
361 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
364 for (
unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
370 jacobian[colIdx][colIdx] = 0;
371 for (
const auto& dataJ : connectivityMap[colIdx])
372 jacobian[dataJ.globalJ][colIdx] = 0;
383 {
return elementColor_[idx]; }
386 {
return elementColor_[idx]; }
390 std::vector<EntityColor> elementColor_;
397template<
class Assembler>
403 using ParentType::ParentType;
409 template<
class Engine>
418template<
class Assembler>
421 using Scalar =
typename Assembler::Scalar;
422 using GridGeometry =
typename Assembler::GridGeometry;
423 using JacobianMatrix =
typename Assembler::JacobianMatrix;
424 using VertexMapper =
typename GridGeometry::VertexMapper;
426 using DiscretizationMethod =
typename GridGeometry::DiscretizationMethod;
439 const auto& gridGeometry = assembler.gridGeometry();
440 totalElems_ = gridGeometry.elementMapper().size();
441 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
460 const std::vector<Scalar>& distanceFromLastLinearization,
463 greenElems_ = engine_.
computeColors(assembler, distanceFromLastLinearization, threshold);
479 template <
class Communication>
480 void report(
const Communication& comm, std::ostream& outStream)
483 greenElems_ = comm.sum(greenElems_);
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_);
498 typename std::enable_if_t<enable, int> = 0>
500 {
return engine_.vertexColor(idx); }
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 resetColors()
Definition: partialreassembler.hh:275
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
void resetColors()
Definition: partialreassembler.hh:377
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.
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
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.