version 3.11-dev
assembly/fvassembler.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// SPDX-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FV_ASSEMBLER_HH
13#define DUMUX_FV_ASSEMBLER_HH
14
15#include <vector>
16#include <deque>
17#include <type_traits>
18#include <memory>
19
20#include <dune/istl/matrixindexset.hh>
21
27
31
35
38
39#include "cvfelocalassembler.hh"
40#include "cclocalassembler.hh"
41#include "fclocalassembler.hh"
42
43namespace Dumux::Detail {
44
45template<class DiscretizationMethod>
47
48template<class DM>
49struct LocalAssemblerChooser<DiscretizationMethods::CVFE<DM>>
50{
51 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
53};
54
55template<>
56struct LocalAssemblerChooser<DiscretizationMethods::CCMpfa>
57{
58 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
60};
61
62template<>
63struct LocalAssemblerChooser<DiscretizationMethods::CCTpfa>
64{
65 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
67};
68
69template<>
70struct LocalAssemblerChooser<DiscretizationMethods::FCStaggered>
71{
72 template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
74};
75
76template<class TypeTag, class Impl, DiffMethod diffMethod, bool isImplicit>
79>::template type<TypeTag, Impl, diffMethod, isImplicit>;
80
82template<class P>
83using ProblemConstraintsDetector = decltype(std::declval<P>().constraints());
84
85template<class P>
86constexpr inline bool hasGlobalConstraints()
87{ return Dune::Std::is_detected<ProblemConstraintsDetector, P>::value; }
88
89} // end namespace Dumux::Detail
90
91namespace Dumux {
92
101template<class TypeTag, DiffMethod diffMethod, bool isImplicit = true, class LocalResidual = GetPropType<TypeTag, Properties::LocalResidual>>
103{
105 using GridView = typename GridGeo::GridView;
106 using Element = typename GridView::template Codim<0>::Entity;
107 using ElementSeed = typename GridView::Grid::template Codim<0>::EntitySeed;
109
110 static constexpr bool isBox = GridGeo::discMethod == DiscretizationMethods::box;
111
114
115public:
120
122
123 using GridGeometry = GridGeo;
125
131 FVAssembler(std::shared_ptr<const Problem> problem,
132 std::shared_ptr<const GridGeometry> gridGeometry,
133 std::shared_ptr<GridVariables> gridVariables)
134 : problem_(problem)
135 , gridGeometry_(gridGeometry)
136 , gridVariables_(gridVariables)
137 , timeLoop_()
138 , isStationaryProblem_(true)
139 {
140 static_assert(isImplicit, "Explicit assembler for stationary problem doesn't make sense!");
142 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
144 && getParam<bool>("Assembly.Multithreading", true);
145
146 maybeComputeColors_();
147 }
148
154 FVAssembler(std::shared_ptr<const Problem> problem,
155 std::shared_ptr<const GridGeometry> gridGeometry,
156 std::shared_ptr<GridVariables> gridVariables,
157 std::shared_ptr<const TimeLoop> timeLoop,
158 const SolutionVector& prevSol)
159 : problem_(problem)
160 , gridGeometry_(gridGeometry)
161 , gridVariables_(gridVariables)
162 , timeLoop_(timeLoop)
163 , prevSol_(&prevSol)
164 , isStationaryProblem_(!timeLoop)
165 {
167 && Grid::Capabilities::supportsMultithreading(gridGeometry_->gridView())
169 && getParam<bool>("Assembly.Multithreading", true);
170
171 maybeComputeColors_();
172 }
173
178 template<class PartialReassembler = DefaultPartialReassembler>
179 void assembleJacobianAndResidual(const SolutionVector& curSol, const PartialReassembler* partialReassembler = nullptr)
180 {
181 checkAssemblerState_();
182 resetJacobian_(partialReassembler);
183 resetResidual_();
184
185 assemble_([&](const Element& element)
186 {
187 LocalAssembler localAssembler(*this, element, curSol);
188 localAssembler.assembleJacobianAndResidual(*jacobian_, *residual_, *gridVariables_, partialReassembler);
189 });
190
191 enforcePeriodicConstraints_(*jacobian_, *residual_, curSol, *gridGeometry_);
192
193 auto applyDirichletConstraint = [&] (const auto& dofIdx,
194 const auto& values,
195 const auto eqIdx,
196 const auto pvIdx)
197 {
198 (*residual_)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
199
200 auto& row = (*jacobian_)[dofIdx];
201 for (auto col = row.begin(); col != row.end(); ++col)
202 row[col.index()][eqIdx] = 0.0;
203
204 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
205 };
206 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
207 }
208
213 {
214 checkAssemblerState_();
215 resetJacobian_();
216
217 assemble_([&](const Element& element)
218 {
219 LocalAssembler localAssembler(*this, element, curSol);
220 localAssembler.assembleJacobian(*jacobian_, *gridVariables_);
221 });
222
223 enforcePeriodicConstraints_(*jacobian_, curSol, *gridGeometry_);
224
225 auto applyDirichletConstraint = [&] (const auto& dofIdx,
226 const auto& values,
227 const auto eqIdx,
228 const auto pvIdx)
229 {
230 auto& row = (*jacobian_)[dofIdx];
231 for (auto col = row.begin(); col != row.end(); ++col)
232 row[col.index()][eqIdx] = 0.0;
233
234 (*jacobian_)[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
235 };
236 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
237 }
238
241 {
242 resetResidual_();
243 assembleResidual(*residual_, curSol);
244 }
245
247 void assembleResidual(ResidualType& r, const SolutionVector& curSol) const
248 {
249 checkAssemblerState_();
250
251 assemble_([&](const Element& element)
252 {
253 LocalAssembler localAssembler(*this, element, curSol);
254 localAssembler.assembleResidual(r);
255 });
256
257 enforcePeriodicConstraints_(*residual_, curSol, *gridGeometry_);
258
259 auto applyDirichletConstraint = [&] (const auto& dofIdx,
260 const auto& values,
261 const auto eqIdx,
262 const auto pvIdx)
263 {
264 r[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
265 };
266 enforceProblemConstraints_(*problem_, *gridGeometry_, applyDirichletConstraint);
267 }
268
274 void setLinearSystem(std::shared_ptr<JacobianMatrix> A,
275 std::shared_ptr<ResidualType> r)
276 {
277 jacobian_ = A;
278 residual_ = r;
279
280 // check and/or set the BCRS matrix's build mode
281 if (jacobian_->buildMode() == JacobianMatrix::BuildMode::unknown)
282 jacobian_->setBuildMode(JacobianMatrix::random);
283 else if (jacobian_->buildMode() != JacobianMatrix::BuildMode::random)
284 DUNE_THROW(Dune::NotImplemented, "Only BCRS matrices with random build mode are supported at the moment");
285
286 setResidualSize_();
287 setJacobianPattern_();
288 }
289
295 {
296 jacobian_ = std::make_shared<JacobianMatrix>();
297 jacobian_->setBuildMode(JacobianMatrix::random);
298 residual_ = std::make_shared<ResidualType>();
299
300 setResidualSize_();
301 setJacobianPattern_();
302 }
303
308 {
309 setResidualSize_();
310 setJacobianPattern_();
311 maybeComputeColors_();
312 }
313
315 std::size_t numDofs() const
316 { return gridGeometry_->numDofs(); }
317
319 const Problem& problem() const
320 { return *problem_; }
321
324 { return *gridGeometry_; }
325
327 const GridView& gridView() const
328 { return gridGeometry().gridView(); }
329
332 { return *gridVariables_; }
333
336 { return *gridVariables_; }
337
340 { return *jacobian_; }
341
344 { return *residual_; }
345
347 const SolutionVector& prevSol() const
348 { return *prevSol_; }
349
354 void setTimeLoop(std::shared_ptr<const TimeLoop> timeLoop)
355 { timeLoop_ = timeLoop; isStationaryProblem_ = !static_cast<bool>(timeLoop); }
356
362 { prevSol_ = &u; }
363
368 { return isStationaryProblem_; }
369
373 LocalResidual localResidual() const
374 { return LocalResidual(problem_.get(), timeLoop_.get()); }
375
380 {
381 gridVariables().update(cursol);
382 }
383
387 void resetTimeStep(const SolutionVector &cursol)
388 {
389 gridVariables().resetTimeStep(cursol);
390 }
391
392private:
396 void setJacobianPattern_()
397 {
398 // resize the jacobian and the residual
399 const auto numDofs = this->numDofs();
400 jacobian_->setSize(numDofs, numDofs);
401
402 // create occupation pattern of the jacobian
403 const auto occupationPattern = getJacobianPattern<isImplicit>(gridGeometry());
404
405 // export pattern to jacobian
406 occupationPattern.exportIdx(*jacobian_);
407 }
408
410 void setResidualSize_()
411 { residual_->resize(numDofs()); }
412
414 void maybeComputeColors_()
415 {
416 if (enableMultithreading_)
417 elementSets_ = computeColoring(gridGeometry()).sets;
418 }
419
420 // reset the residual vector to 0.0
421 void resetResidual_()
422 {
423 if(!residual_)
424 {
425 residual_ = std::make_shared<ResidualType>();
426 setResidualSize_();
427 }
428
429 (*residual_) = 0.0;
430 }
431
432 // reset the Jacobian matrix to 0.0
433 template <class PartialReassembler = DefaultPartialReassembler>
434 void resetJacobian_(const PartialReassembler *partialReassembler = nullptr)
435 {
436 if(!jacobian_)
437 {
438 jacobian_ = std::make_shared<JacobianMatrix>();
439 jacobian_->setBuildMode(JacobianMatrix::random);
440 setJacobianPattern_();
441 }
442
443 if (partialReassembler)
444 partialReassembler->resetJacobian(*this);
445 else
446 *jacobian_ = 0.0;
447 }
448
449 // check if the assembler is in a correct state for assembly
450 void checkAssemblerState_() const
451 {
452 if (!isStationaryProblem_ && !prevSol_)
453 DUNE_THROW(Dune::InvalidStateException, "Assembling instationary problem but previous solution was not set!");
454 }
455
461 template<typename AssembleElementFunc>
462 void assemble_(AssembleElementFunc&& assembleElement) const
463 {
464 // a state that will be checked on all processes
465 bool succeeded = false;
466
467 // try assembling using the local assembly function
468 try
469 {
470 if (enableMultithreading_)
471 {
472 assert(elementSets_.size() > 0);
473
474 // make this element loop run in parallel
475 // for this we have to color the elements so that we don't get
476 // race conditions when writing into the global matrix
477 // each color can be assembled using multiple threads
478 for (const auto& elements : elementSets_)
479 {
480 Dumux::parallelFor(elements.size(), [&](const std::size_t i)
481 {
482 const auto element = gridView().grid().entity(elements[i]);
483 assembleElement(element);
484 });
485 }
486 }
487 else
488 for (const auto& element : elements(gridView()))
489 assembleElement(element);
490
491 // if we get here, everything worked well on this process
492 succeeded = true;
493 }
494 // throw exception if a problem occurred
495 catch (NumericalProblem &e)
496 {
497 std::cout << "rank " << gridView().comm().rank()
498 << " caught an exception while assembling:" << e.what()
499 << "\n";
500 succeeded = false;
501 }
502
503 // make sure everything worked well on all processes
504 if (gridView().comm().size() > 1)
505 succeeded = gridView().comm().min(succeeded);
506
507 // if not succeeded rethrow the error on all processes
508 if (!succeeded)
509 DUNE_THROW(NumericalProblem, "A process did not succeed in linearizing the system");
510 }
511
512 template<class GG>
513 void enforcePeriodicConstraints_(JacobianMatrix& jac, ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry) const
514 {
515 if constexpr (Detail::hasPeriodicDofMap<GG>())
516 {
517 for (const auto& m : gridGeometry.periodicDofMap())
518 {
519 if (m.first < m.second)
520 {
521 // add the second row to the first
522 res[m.first] += res[m.second];
523 const auto end = jac[m.second].end();
524 for (auto it = jac[m.second].begin(); it != end; ++it)
525 jac[m.first][it.index()] += (*it);
526
527 // enforce constraint in second row
528 res[m.second] = curSol[m.second] - curSol[m.first];
529
530 // set derivatives accordingly in jacobian, i.e. id for m.second and -id for m.first
531 auto setMatrixBlock = [] (auto& matrixBlock, double diagValue)
532 {
533 for (int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
534 matrixBlock[eIdx][eIdx] = diagValue;
535 };
536
537 for (auto it = jac[m.second].begin(); it != end; ++it)
538 {
539 auto& matrixBlock = *it;
540 matrixBlock = 0.0;
541
542 assert(matrixBlock.N() == matrixBlock.M());
543 if(it.index() == m.second)
544 setMatrixBlock(matrixBlock, 1.0);
545
546 if(it.index() == m.first)
547 setMatrixBlock(matrixBlock, -1.0);
548
549 }
550 }
551 }
552 }
553 }
554
555 template<class GG>
556 void enforcePeriodicConstraints_(JacobianMatrix& jac, const SolutionVector& curSol, const GG& gridGeometry) const
557 {
558 if constexpr (Detail::hasPeriodicDofMap<GG>())
559 {
560 for (const auto& m : gridGeometry.periodicDofMap())
561 {
562 if (m.first < m.second)
563 {
564 // add the second row to the first
565 const auto end = jac[m.second].end();
566 for (auto it = jac[m.second].begin(); it != end; ++it)
567 jac[m.first][it.index()] += (*it);
568
569 // set derivatives accordingly in jacobian, i.e. id for m.second and -id for m.first
570 auto setMatrixBlock = [] (auto& matrixBlock, double diagValue)
571 {
572 for (int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
573 matrixBlock[eIdx][eIdx] = diagValue;
574 };
575
576 for (auto it = jac[m.second].begin(); it != end; ++it)
577 {
578 auto& matrixBlock = *it;
579 matrixBlock = 0.0;
580
581 assert(matrixBlock.N() == matrixBlock.M());
582 if(it.index() == m.second)
583 setMatrixBlock(matrixBlock, 1.0);
584
585 if(it.index() == m.first)
586 setMatrixBlock(matrixBlock, -1.0);
587
588 }
589 }
590 }
591 }
592 }
593
594 template<class GG>
595 void enforcePeriodicConstraints_(ResidualType& res, const SolutionVector& curSol, const GG& gridGeometry) const
596 {
597 if constexpr (Detail::hasPeriodicDofMap<GG>())
598 {
599 for (const auto& m : gridGeometry.periodicDofMap())
600 {
601 if (m.first < m.second)
602 {
603 // add the second row to the first
604 res[m.first] += res[m.second];
605
606 // enforce constraint in second row
607 res[m.second] = curSol[m.second] - curSol[m.first];
608 }
609 }
610 }
611 }
612
613 template<class Problem, class GG, typename ApplyFunction>
614 void enforceProblemConstraints_(const Problem& problem, const GG& gridGeometry, const ApplyFunction& applyDirichletConstraint) const
615 {
616 if constexpr (Detail::hasGlobalConstraints<Problem>())
617 {
618 for (const auto& constraintData : problem.constraints())
619 {
620 const auto& constraintInfo = constraintData.constraintInfo();
621 const auto& values = constraintData.values();
622 const auto dofIdx = constraintData.dofIndex();
623 // set the constraint in residual and jacobian
624 for (int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
625 {
626 if (constraintInfo.isConstraintEquation(eqIdx))
627 {
628 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
629 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
630 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
631 }
632 }
633 }
634 }
635 }
636
638 std::shared_ptr<const Problem> problem_;
639
641 std::shared_ptr<const GridGeometry> gridGeometry_;
642
644 std::shared_ptr<GridVariables> gridVariables_;
645
647 std::shared_ptr<const TimeLoop> timeLoop_;
648
650 const SolutionVector* prevSol_ = nullptr;
651
653 bool isStationaryProblem_;
654
656 std::shared_ptr<JacobianMatrix> jacobian_;
657 std::shared_ptr<ResidualType> residual_;
658
660 bool enableMultithreading_ = false;
661 std::deque<std::vector<ElementSeed>> elementSets_;
662};
663
664} // namespace Dumux
665
666#endif
An assembler for Jacobian and residual contribution per element (cell-centered methods)
Definition: assembly/cclocalassembler.hh:124
An assembler for Jacobian and residual contribution per element (CVFE methods)
Definition: assembly/cvfelocalassembler.hh:317
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: assembly/fvassembler.hh:103
void setPreviousSolution(const SolutionVector &u)
Sets the solution from which to start the time integration. Has to be called prior to assembly for ti...
Definition: assembly/fvassembler.hh:361
GridGeo GridGeometry
Definition: assembly/fvassembler.hh:123
GetPropType< TypeTag, Properties::GridVariables > GridVariables
Definition: assembly/fvassembler.hh:121
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: assembly/fvassembler.hh:154
GridVariables & gridVariables()
The global grid variables.
Definition: assembly/fvassembler.hh:331
std::size_t numDofs() const
Returns the number of degrees of freedom.
Definition: assembly/fvassembler.hh:315
const GridVariables & gridVariables() const
The global grid variables.
Definition: assembly/fvassembler.hh:335
LocalResidual localResidual() const
Create a local residual object (used by the local assembler)
Definition: assembly/fvassembler.hh:373
GetPropType< TypeTag, Properties::SolutionVector > SolutionVector
Definition: assembly/fvassembler.hh:118
void assembleJacobianAndResidual(const SolutionVector &curSol, const PartialReassembler *partialReassembler=nullptr)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: assembly/fvassembler.hh:179
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: assembly/fvassembler.hh:307
const Problem & problem() const
The problem.
Definition: assembly/fvassembler.hh:319
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: assembly/fvassembler.hh:367
const GridGeometry & gridGeometry() const
The global finite volume geometry.
Definition: assembly/fvassembler.hh:323
void assembleResidual(ResidualType &r, const SolutionVector &curSol) const
assemble a residual r
Definition: assembly/fvassembler.hh:247
typename Detail::NativeDuneVectorType< SolutionVector >::type ResidualType
Definition: assembly/fvassembler.hh:119
void setTimeLoop(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: assembly/fvassembler.hh:354
void resetTimeStep(const SolutionVector &cursol)
Reset the gridVariables.
Definition: assembly/fvassembler.hh:387
void assembleJacobian(const SolutionVector &curSol)
Assembles only the global Jacobian of the residual.
Definition: assembly/fvassembler.hh:212
void updateGridVariables(const SolutionVector &cursol)
Update the grid variables.
Definition: assembly/fvassembler.hh:379
JacobianMatrix & jacobian()
The jacobian matrix.
Definition: assembly/fvassembler.hh:339
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< ResidualType > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: assembly/fvassembler.hh:274
GetPropType< TypeTag, Properties::JacobianMatrix > JacobianMatrix
Definition: assembly/fvassembler.hh:117
GetPropType< TypeTag, Properties::Problem > Problem
Definition: assembly/fvassembler.hh:124
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: assembly/fvassembler.hh:294
const SolutionVector & prevSol() const
The solution of the previous time step.
Definition: assembly/fvassembler.hh:347
FVAssembler(std::shared_ptr< const Problem > problem, std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< GridVariables > gridVariables)
The constructor for stationary problems.
Definition: assembly/fvassembler.hh:131
const GridView & gridView() const
The gridview.
Definition: assembly/fvassembler.hh:327
ResidualType & residual()
The residual vector (rhs)
Definition: assembly/fvassembler.hh:343
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: assembly/fvassembler.hh:240
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: assembly/fvassembler.hh:116
An assembler for Jacobian and residual contribution per element (Face-centered methods)
Definition: fclocalassembler.hh:288
detects which entries in the Jacobian have to be recomputed
Definition: partialreassembler.hh:420
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:84
Coloring schemes for shared-memory-parallel assembly.
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper to extract native Dune vector types from particular Dumux types.
An assembler for Jacobian and residual contribution per element (cell-centered methods)
An assembler for Jacobian and residual contribution per element (CVFE methods)
An assembler for Jacobian and residual contribution per element (face-centered staggered methods)
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
void parallelFor(const std::size_t count, const FunctorType &functor)
A parallel for loop (multithreading)
Definition: parallel_for.hh:160
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
Multithreading in Dumux.
Definition: cvfelocalresidual.hh:30
constexpr bool hasGlobalConstraints()
Definition: assembly/fvassembler.hh:86
decltype(std::declval< P >().constraints()) ProblemConstraintsDetector
helper struct detecting if problem has a constraints() function
Definition: assembly/fvassembler.hh:83
typename LocalAssemblerChooser< typename GetPropType< TypeTag, Properties::GridGeometry >::DiscretizationMethod >::template type< TypeTag, Impl, diffMethod, isImplicit > LocalAssemblerChooser_t
Definition: assembly/fvassembler.hh:79
constexpr Box box
Definition: method.hh:156
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:34
Definition: adapt.hh:17
auto computeColoring(const GridGeometry &gg, int verbosity=1)
Compute iterable lists of element seeds partitioned by color.
Definition: coloring.hh:241
Parallel for loop (multithreading)
Provides a helper class for nonoverlapping decomposition.
Type traits to detect periodicity support.
Definition: assembly/fvassembler.hh:46
typename NativeDuneVectorTypeImpl< V, Dune::Std::is_detected< Detail::DuneVectors::StateDetector, V >{} >::type type
Definition: dunevectors.hh:59
Traits specifying if a given discretization tag supports coloring.
Definition: coloring.hh:294
Type traits to be used with vector types.