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>
44 template<
typename... Args>
46 { DUNE_THROW(Dune::InvalidStateException,
"DefaultPartialReassembler should be never constructed!"); }
48 template<
typename... Args>
51 template<
typename... Args>
54 template<
typename... Args>
57 template<
typename... Args>
71template<
class Assembler, DiscretizationMethod discMethod>
76 { DUNE_THROW(Dune::NotImplemented,
"PartialReassembler for this discretization method!"); }
84 template<
typename... Args>
87 template<
typename... Args>
90 template<
typename... Args>
98template<
class Assembler>
101 using Scalar =
typename Assembler::Scalar;
102 using GridGeometry =
typename Assembler::GridGeometry;
103 using JacobianMatrix =
typename Assembler::JacobianMatrix;
104 using VertexMapper =
typename GridGeometry::VertexMapper;
105 static constexpr int dim = GridGeometry::GridView::dimension;
109 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
110 , vertexColor_(assembler.gridGeometry().vertexMapper().size(),
EntityColor::
red)
115 const std::vector<Scalar>& distanceFromLastLinearization,
118 const auto& gridGeometry = assembler.gridGeometry();
119 const auto& gridView = gridGeometry.gridView();
120 const auto& elementMapper = gridGeometry.elementMapper();
121 const auto& vertexMapper = gridGeometry.vertexMapper();
127 for (
unsigned int i = 0; i < vertexColor_.size(); ++i)
130 if (distanceFromLastLinearization[i] > threshold)
137 for (
const auto& element : elements(gridView))
142 int numVertices = element.subEntities(dim);
144 for (
int i = 0; i < numVertices; ++i) {
145 int globalI = vertexMapper.subIndex(element, i, dim);
153 int eIdx = elementMapper.index(element);
162 for (
const auto& element : elements(gridView))
164 int eIdx = elementMapper.index(element);
169 int numVertices = element.subEntities(dim);
171 for (
int i = 0; i < numVertices; ++i) {
172 int globalI = vertexMapper.subIndex(element, i, dim);
185 minHandle(vertexColor_, vertexMapper);
186 gridView.communicate(minHandle,
187 Dune::InteriorBorder_InteriorBorder_Interface,
188 Dune::ForwardCommunication);
191 for (
const auto& element : elements(gridView))
193 int eIdx = elementMapper.index(element);
199 bool isOrange =
false;
200 int numVertices = element.subEntities(dim);
202 for (
int i = 0; i < numVertices; ++i) {
203 int globalI = vertexMapper.subIndex(element, i, dim);
218 for (
const auto& element : elements(gridView))
220 int eIdx = elementMapper.index(element);
225 int numVertices = element.subEntities(dim);
227 for (
int i = 0; i < numVertices; ++i) {
228 int globalI = vertexMapper.subIndex(element, i, dim);
239 maxHandle(vertexColor_, vertexMapper);
240 gridView.communicate(maxHandle,
241 Dune::InteriorBorder_InteriorBorder_Interface,
242 Dune::ForwardCommunication);
245 for (
unsigned int i=0; i < vertexColor_.size(); ++i) {
255 return std::count_if(elementColor_.begin(), elementColor_.end(),
256 [](
EntityColor c){ return c == EntityColor::green; });
261 auto& jacobian = assembler.jacobian();
264 for (
unsigned int rowIdx = 0; rowIdx < jacobian.N(); ++rowIdx)
270 auto colIt = jacobian[rowIdx].begin();
271 const auto& colEndIt = jacobian[rowIdx].end();
272 for (; colIt != colEndIt; ++colIt) {
286 {
return elementColor_[idx]; }
289 {
return vertexColor_[idx]; }
292 {
return vertexColor_[idx]; }
296 std::vector<EntityColor> elementColor_;
297 std::vector<EntityColor> vertexColor_;
304template<
class Assembler>
307 using Scalar =
typename Assembler::Scalar;
308 using GridGeometry =
typename Assembler::GridGeometry;
309 using JacobianMatrix =
typename Assembler::JacobianMatrix;
313 : elementColor_(assembler.gridGeometry().elementMapper().size(),
EntityColor::
red)
318 const std::vector<Scalar>& distanceFromLastLinearization,
321 const auto& gridGeometry = assembler.gridGeometry();
322 const auto& gridView = gridGeometry.gridView();
323 const auto& elementMapper = gridGeometry.elementMapper();
326 for (
const auto& element : elements(gridView))
328 int eIdx = elementMapper.index(element);
330 if (distanceFromLastLinearization[eIdx] > threshold)
343 const auto& connectivityMap = gridGeometry.connectivityMap();
344 for (
unsigned eIdx = 0; eIdx < elementColor_.size(); ++eIdx)
349 if (distanceFromLastLinearization[eIdx] > threshold)
351 for (
const auto& connectedDof : connectivityMap[eIdx])
357 return std::count_if(elementColor_.begin(), elementColor_.end(),
358 [](
EntityColor c){return c == EntityColor::green;});
364 auto& jacobian = assembler.jacobian();
365 const auto& connectivityMap = assembler.gridGeometry().connectivityMap();
368 for (
unsigned int colIdx = 0; colIdx < jacobian.M(); ++colIdx)
374 jacobian[colIdx][colIdx] = 0;
375 for (
const auto& dataJ : connectivityMap[colIdx])
376 jacobian[dataJ.globalJ][colIdx] = 0;
387 {
return elementColor_[idx]; }
390 {
return elementColor_[idx]; }
394 std::vector<EntityColor> elementColor_;
401template<
class Assembler>
407 using ParentType::ParentType;
413 template<
class Engine>
422template<
class Assembler>
425 using Scalar =
typename Assembler::Scalar;
426 using GridGeometry =
typename Assembler::GridGeometry;
427 using JacobianMatrix =
typename Assembler::JacobianMatrix;
428 using VertexMapper =
typename GridGeometry::VertexMapper;
443 const auto& gridGeometry = assembler.gridGeometry();
444 totalElems_ = gridGeometry.elementMapper().size();
445 totalElems_ = gridGeometry.gridView().comm().sum(totalElems_);
464 const std::vector<Scalar>& distanceFromLastLinearization,
467 greenElems_ = engine_.computeColors(assembler, distanceFromLastLinearization, threshold);
472 engine_.resetColors();
477 engine_.resetJacobian(assembler);
483 template <
class Communication>
484 void report(
const Communication& comm, std::ostream& outStream)
487 greenElems_ = comm.sum(greenElems_);
489 auto reassembledElems = totalElems_ - greenElems_;
490 auto width = std::to_string(totalElems_).size();
491 outStream <<
", reassembled " << std::setw(width)
492 << reassembledElems <<
" (" << std::setw(3)
493 << 100*reassembledElems/totalElems_ <<
"%) elements";
497 {
return engine_.elementColor(idx); }
500 {
return engine_.dofColor(idx); }
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.
Provides data handles for parallel communication which operate on vertices.
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ box
Definition method.hh:38
@ ccmpfa
Definition method.hh:38
@ cctpfa
Definition method.hh:38
EntityColor
The colors of elements and vertices required for partial Jacobian reassembly.
Definition entitycolor.hh:35
@ green
does not need to be reassembled
Definition entitycolor.hh:52
@ red
distance from last linearization is above the tolerance
Definition entitycolor.hh:37
@ yellow
neighboring entity is red
Definition entitycolor.hh:40
@ orange
Definition entitycolor.hh:49
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
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
EntityColor dofColor(size_t idx) const
Definition partialreassembler.hh:60
void resetJacobian(Args &&... args) const
Definition partialreassembler.hh:52
void computeColors(Args &&... args)
Definition partialreassembler.hh:55
EntityColor vertexColor(size_t idx) const
Definition partialreassembler.hh:66
void report(Args &&... args)
Definition partialreassembler.hh:49
DefaultPartialReassembler(Args &&... args)
Definition partialreassembler.hh:45
void resetColors(Args &&... args)
Definition partialreassembler.hh:58
EntityColor elementColor(size_t idx) const
Definition partialreassembler.hh:63
the partial reassembler engine specialized for discretization methods
Definition partialreassembler.hh:73
std::size_t computeColors(Args &&... args)
Definition partialreassembler.hh:85
EntityColor elementColor(size_t idx) const
Definition partialreassembler.hh:78
EntityColor dofColor(size_t idx) const
Definition partialreassembler.hh:81
void resetColors(Args &&... args)
Definition partialreassembler.hh:91
PartialReassemblerEngine(const Assembler &)
Definition partialreassembler.hh:75
void resetJacobian(Args &&... args) const
Definition partialreassembler.hh:88
EntityColor dofColor(size_t idx) const
Definition partialreassembler.hh:291
void resetColors()
Definition partialreassembler.hh:279
PartialReassemblerEngine(const Assembler &assembler)
Definition partialreassembler.hh:108
EntityColor vertexColor(size_t idx) const
Definition partialreassembler.hh:288
EntityColor elementColor(size_t idx) const
Definition partialreassembler.hh:285
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition partialreassembler.hh:114
void resetJacobian(Assembler &assembler) const
Definition partialreassembler.hh:259
void resetJacobian(Assembler &assembler) const
Definition partialreassembler.hh:362
void resetColors()
Definition partialreassembler.hh:381
EntityColor elementColor(size_t idx) const
Definition partialreassembler.hh:386
PartialReassemblerEngine(const Assembler &assembler)
Definition partialreassembler.hh:312
std::size_t computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Definition partialreassembler.hh:317
EntityColor dofColor(size_t idx) const
Definition partialreassembler.hh:389
helper struct to determine whether the an engine class has vertex colors
Definition partialreassembler.hh:412
auto operator()(Engine &&e) -> decltype(e.vertexColor(0))
Definition partialreassembler.hh:414
PartialReassembler(const Assembler &assembler)
constructor
Definition partialreassembler.hh:439
EntityColor vertexColor(size_t idx) const
Definition partialreassembler.hh:504
void resetJacobian(Assembler &assembler) const
Definition partialreassembler.hh:475
void resetColors()
Definition partialreassembler.hh:470
void computeColors(const Assembler &assembler, const std::vector< Scalar > &distanceFromLastLinearization, Scalar threshold)
Determine the colors of entities for partial reassembly.
Definition partialreassembler.hh:463
void report(const Communication &comm, std::ostream &outStream)
called by the assembler after successful assembly
Definition partialreassembler.hh:484
EntityColor dofColor(size_t idx) const
Definition partialreassembler.hh:499
EntityColor elementColor(size_t idx) const
Definition partialreassembler.hh:496
Data handle for parallel communication which takes the maximum of all values that are attached to ver...
Definition vertexhandles.hh:100
Provides data handle for parallel communication which takes the minimum of all values that are attach...
Definition vertexhandles.hh:160