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>
46 template<
typename... Args>
48 { DUNE_THROW(Dune::InvalidStateException,
"DefaultPartialReassembler should be never constructed!"); }
50 template<
typename... Args>
53 template<
typename... Args>
56 template<
typename... Args>
59 template<
typename... Args>
73template<
class Assembler,
class DiscretizationMethod>
78 { DUNE_THROW(Dune::NotImplemented,
"PartialReassembler for this discretization method!"); }
86 template<
typename... Args>
89 template<
typename... Args>
92 template<
typename... Args>
100template<
class Assembler>
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;
111 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
112 , vertexColor_(assembler.gridGeometry().vertexMapper().size(),
EntityColor::
red)
117 const std::vector<Scalar>& distanceFromLastLinearization,
120 const auto& gridGeometry = assembler.gridGeometry();
121 const auto& gridView = gridGeometry.gridView();
122 const auto& elementMapper = gridGeometry.elementMapper();
123 const auto& vertexMapper = gridGeometry.vertexMapper();
129 for (
unsigned int i = 0; i < vertexColor_.size(); ++i)
132 if (distanceFromLastLinearization[i] > threshold)
139 for (
const auto& element : elements(gridView))
144 int numVertices = element.subEntities(dim);
146 for (
int i = 0; i < numVertices; ++i) {
147 int globalI = vertexMapper.subIndex(element, i, dim);
155 int eIdx = elementMapper.index(element);
164 for (
const auto& element : elements(gridView))
166 int eIdx = elementMapper.index(element);
171 int numVertices = element.subEntities(dim);
173 for (
int i = 0; i < numVertices; ++i) {
174 int globalI = vertexMapper.subIndex(element, i, dim);
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);
193 DUNE_THROW(Dune::InvalidStateException,
"Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim <<
"-entities.");
196 for (
const auto& element : elements(gridView))
198 int eIdx = elementMapper.index(element);
204 bool isOrange =
false;
205 int numVertices = element.subEntities(dim);
207 for (
int i = 0; i < numVertices; ++i) {
208 int globalI = vertexMapper.subIndex(element, i, dim);
223 for (
const auto& element : elements(gridView))
225 int eIdx = elementMapper.index(element);
230 int numVertices = element.subEntities(dim);
232 for (
int i = 0; i < numVertices; ++i) {
233 int globalI = vertexMapper.subIndex(element, i, dim);
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);
250 DUNE_THROW(Dune::InvalidStateException,
"Cannot call computeColors on multiple processes for a grid that cannot communicate codim-" << dim <<
"-entities.");
253 for (
unsigned int i=0; i < vertexColor_.size(); ++i) {
263 return std::count_if(elementColor_.begin(), elementColor_.end(),
264 [](
EntityColor c){ return c == EntityColor::green; });
269 auto& jacobian = assembler.jacobian();
272 for (
unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
278 auto colIt = jacobian[rowIdx].begin();
279 const auto& colEndIt = jacobian[rowIdx].end();
280 for (; colIt != colEndIt; ++colIt) {
294 {
return elementColor_[idx]; }
297 {
return vertexColor_[idx]; }
300 {
return vertexColor_[idx]; }
304 std::vector<EntityColor> elementColor_;
305 std::vector<EntityColor> vertexColor_;
312template<
class Assembler>
315 using Scalar =
typename Assembler::Scalar;
316 using GridGeometry =
typename Assembler::GridGeometry;
317 using JacobianMatrix =
typename Assembler::JacobianMatrix;
321 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
326 const std::vector<Scalar>& distanceFromLastLinearization,
329 const auto& gridGeometry = assembler.gridGeometry();
330 const auto& gridView = gridGeometry.gridView();
331 const auto& elementMapper = gridGeometry.elementMapper();
334 for (
const auto& element : elements(gridView))
336 int eIdx = elementMapper.index(element);
338 if (distanceFromLastLinearization[eIdx] > threshold)
351 const auto& connectivityMap = gridGeometry.connectivityMap();
352 for (
unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
357 if (distanceFromLastLinearization[eIdx] > threshold)
359 for (
const auto& connectedDof : connectivityMap[eIdx])
365 return std::count_if(elementColor_.begin(), elementColor_.end(),
366 [](
EntityColor c){return c == EntityColor::green;});
372 auto& jacobian = assembler.jacobian();
373 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
376 for (
unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
382 jacobian[colIdx][colIdx] = 0;
383 for (
const auto& dataJ : connectivityMap[colIdx])
384 jacobian[dataJ.globalJ][colIdx] = 0;
395 {
return elementColor_[idx]; }
398 {
return elementColor_[idx]; }
402 std::vector<EntityColor> elementColor_;
409template<
class Assembler>
415 using ParentType::ParentType;
421 template<
class Engine>
430template<
class Assembler>
433 using Scalar =
typename Assembler::Scalar;
434 using GridGeometry =
typename Assembler::GridGeometry;
435 using JacobianMatrix =
typename Assembler::JacobianMatrix;
436 using VertexMapper =
typename GridGeometry::VertexMapper;
438 using DiscretizationMethod =
typename GridGeometry::DiscretizationMethod;
451 const auto& gridGeometry = assembler.gridGeometry();
452 totalElems_ = gridGeometry.elementMapper().size();
453 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
472 const std::vector<Scalar>& distanceFromLastLinearization,
475 greenElems_ = engine_.
computeColors(assembler, distanceFromLastLinearization, threshold);
491 template <
class Communication>
492 void report(
const Communication& comm, std::ostream& outStream)
495 greenElems_ = comm.sum(greenElems_);
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_);
510 typename std::enable_if_t<enable, int> = 0>
512 {
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.
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 resetColors()
Definition: partialreassembler.hh:287
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
void resetColors()
Definition: partialreassembler.hh:389
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