14#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
15#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
20#include <dune/common/hybridutilities.hh>
21#include <dune/istl/matrixindexset.hh>
42#if HAVE_DUMUX_OLD_STAGGERED
43#include <dumux/multidomain/subdomainstaggeredlocalassembler.hh>
48template<std::
size_t id,
class TypeTag,
class Assembler, DiffMethod DM,
bool implicit>
51namespace Grid::Capabilities {
55template<
class T, std::size_t... I>
79{
return Dune::Std::is_detected<SubProblemConstraintsDetector, P>::value; }
89struct CouplingManagerSupportsMultithreadedAssembly :
public std::false_type
101template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
104 template<std::
size_t id>
105 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
110 using Scalar =
typename MDTraits::Scalar;
112 template<std::
size_t id>
113 using LocalResidual =
typename MDTraits::template SubDomain<id>::LocalResidual;
115 template<std::
size_t id>
116 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
118 template<std::
size_t id>
119 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
121 template<std::
size_t id>
122 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
134 {
return useImplicitAssembly; }
138 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
139 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
140 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
145 template<std::
size_t id>
148 template<
class DiscretizationMethod, std::
size_t id>
149 struct SubDomainAssemblerType;
151 template<std::
size_t id>
152 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
157 template<std::
size_t id>
158 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
160 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
163 template<std::
size_t id,
class DM>
164 struct SubDomainAssemblerType<DiscretizationMethods::CVFE<DM>, id>
166 using type = SubDomainCVFELocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
169 template<std::
size_t id>
170 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
172 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
175 template<std::
size_t id>
176 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
178 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, SubDomainAssemblerView<id>, diffMethod,
isImplicit()>;
181 template<std::
size_t id>
182 using SubDomainAssembler =
typename SubDomainAssemblerType<typename GridGeometry<id>::DiscretizationMethod,
id>::type;
197 , problemTuple_(std::move(
problem))
201 , isStationaryProblem_(true)
202 , warningIssued_(false)
204 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
205 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
210 && getParam<bool>(
"Assembly.Multithreading",
true);
212 maybeComputeColors_();
224 std::shared_ptr<const TimeLoop> timeLoop,
227 , problemTuple_(std::move(
problem))
230 , timeLoop_(timeLoop)
232 , isStationaryProblem_(false)
233 , warningIssued_(false)
235 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
240 && getParam<bool>(
"Assembly.Multithreading",
true);
242 maybeComputeColors_();
251 checkAssemblerState_();
255 using namespace Dune::Hybrid;
256 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
258 auto& jacRow = (*jacobian_)[domainId];
259 auto& subRes = (*residual_)[domainId];
260 const auto& subSol = curSol[domainId];
261 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
263 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
264 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, subSol);
265 enforceProblemConstraints_(domainId, jacRow, subRes, *
gridGeometry, subSol);
281 checkAssemblerState_();
286 using namespace Dune::Hybrid;
287 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
289 auto& subRes = r[domainId];
290 const auto& subSol = curSol[domainId];
291 this->assembleResidual_(domainId, subRes, curSol);
293 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
294 enforcePeriodicConstraints_(domainId, subRes, *
gridGeometry, subSol);
295 enforceProblemConstraints_(domainId, subRes, *
gridGeometry, subSol);
305 std::shared_ptr<ResidualType> r)
311 setJacobianPattern_(*jacobian_);
312 setResidualSize_(*residual_);
321 jacobian_ = std::make_shared<JacobianMatrix>();
322 residual_ = std::make_shared<ResidualType>();
325 setJacobianPattern_(*jacobian_);
326 setResidualSize_(*residual_);
334 using namespace Dune::Hybrid;
335 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
337 forEach(jac[i], [&](
auto& jacBlock)
339 using BlockType = std::decay_t<
decltype(jacBlock)>;
340 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
341 jacBlock.setBuildMode(BlockType::BuildMode::random);
342 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
343 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
353 setJacobianPattern_(*jacobian_);
354 setResidualSize_(*residual_);
355 maybeComputeColors_();
363 using namespace Dune::Hybrid;
364 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
365 { this->
gridVariables(domainId).update(curSol[domainId]); });
373 using namespace Dune::Hybrid;
374 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
375 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
379 template<std::
size_t i>
380 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
381 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
384 template<std::
size_t i>
385 const auto&
problem(Dune::index_constant<i> domainId)
const
386 {
return *std::get<domainId>(problemTuple_); }
389 template<std::
size_t i>
391 {
return *std::get<domainId>(gridGeometryTuple_); }
394 template<std::
size_t i>
395 const auto&
gridView(Dune::index_constant<i> domainId)
const
399 template<std::
size_t i>
401 {
return *std::get<domainId>(gridVariablesTuple_); }
404 template<std::
size_t i>
406 {
return *std::get<domainId>(gridVariablesTuple_); }
414 {
return *jacobian_; }
418 {
return *residual_; }
422 {
return *prevSol_; }
429 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
442 {
return isStationaryProblem_; }
447 template<std::
size_t i>
449 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
461 using namespace Dune::Hybrid;
462 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
464 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
466 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
467 pattern.exportIdx(jac[domainI][domainJ]);
477 using namespace Dune::Hybrid;
478 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
479 { res[domainId].resize(this->
numDofs(domainId)); });
483 void resetResidual_()
487 residual_ = std::make_shared<ResidualType>();
488 setResidualSize_(*residual_);
495 void resetJacobian_()
499 jacobian_ = std::make_shared<JacobianMatrix>();
501 setJacobianPattern_(*jacobian_);
508 void maybeComputeColors_()
510 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
511 if (enableMultithreading_)
516 void checkAssemblerState_()
const
518 if (!isStationaryProblem_ && !prevSol_)
519 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
521 if (isStationaryProblem_ && prevSol_)
522 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
523 <<
" Did you forget to set the timeLoop to make this problem instationary?");
526 template<std::
size_t i,
class JacRow,
class SubRes>
527 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
530 assemble_(domainId, [&](
const auto& element)
532 MultiDomainAssemblerSubDomainView view{*
this, domainId};
533 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
534 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
538 template<std::
size_t i,
class SubRes>
539 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
542 assemble_(domainId, [&](
const auto& element)
544 MultiDomainAssemblerSubDomainView view{*
this, domainId};
545 SubDomainAssembler<i> subDomainAssembler(view, element, curSol, *
couplingManager_);
546 subDomainAssembler.assembleResidual(subRes);
555 template<std::
size_t i,
class AssembleElementFunc>
556 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
559 bool succeeded =
false;
564 if constexpr (CouplingManagerSupportsMultithreadedAssembly<CouplingManager>::value)
566 if (enableMultithreading_)
569 domainId, std::forward<AssembleElementFunc>(assembleElement)
578 for (
const auto& element : elements(
gridView(domainId)))
579 assembleElement(element);
585 catch (NumericalProblem &e)
587 std::cout <<
"rank " <<
gridView(domainId).comm().rank()
588 <<
" caught an exception while assembling:" << e.what()
594 if (
gridView(domainId).comm().size() > 1)
595 succeeded =
gridView(domainId).comm().min(succeeded);
599 DUNE_THROW(NumericalProblem,
"A process did not succeed in linearizing the system");
603 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
604 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
605 Dune::index_constant<j> domainJ)
const
608 auto pattern = getJacobianPattern<isImplicit()>(gg);
614 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
615 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
616 Dune::index_constant<j> domainJ)
const
624 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
625 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
627 if constexpr (Detail::hasPeriodicDofMap<GG>())
631 if (m.first < m.second)
633 auto& jac = jacRow[domainI];
636 res[m.first] += res[m.second];
638 const auto end = jac[m.second].end();
639 for (
auto it = jac[m.second].begin(); it != end; ++it)
640 jac[m.first][it.index()] += (*it);
644 res[m.second] = curSol[m.second] - curSol[m.first];
647 auto setMatrixBlock = [] (
auto& matrixBlock,
double diagValue)
649 for (
int eIdx = 0; eIdx < matrixBlock.N(); ++eIdx)
650 matrixBlock[eIdx][eIdx] = diagValue;
653 for (
auto it = jac[m.second].begin(); it != end; ++it)
655 auto& matrixBlock = *it;
658 assert(matrixBlock.N() == matrixBlock.M());
659 if(it.index() == m.second)
660 setMatrixBlock(matrixBlock, 1.0);
662 if(it.index() == m.first)
663 setMatrixBlock(matrixBlock, -1.0);
667 using namespace Dune::Hybrid;
670 auto& jacCoupling = jacRow[couplingDomainId];
672 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
673 jacCoupling[m.first][it.index()] += (*it);
675 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
684 template<std::
size_t i,
class JacRow,
class Res,
class GG,
class Sol>
685 void enforceProblemConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
687 if constexpr (Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
689 auto& jac = jacRow[domainI];
690 auto applyDirichletConstraint = [&] (
const auto& dofIdx,
695 (res)[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
697 auto& row = jac[dofIdx];
698 for (
auto col = row.begin(); col != row.end(); ++col)
699 row[col.index()][eqIdx] = 0.0;
701 jac[dofIdx][dofIdx][eqIdx][pvIdx] = 1.0;
703 using namespace Dune::Hybrid;
706 auto& jacCoupling = jacRow[couplingDomainId];
708 auto& rowCoupling = jacCoupling[dofIdx];
709 for (
auto colCoupling = rowCoupling.begin(); colCoupling != rowCoupling.end(); ++colCoupling)
710 rowCoupling[colCoupling.index()][eqIdx] = 0.0;
715 for (
const auto& constraintData : this->
problem(domainI).constraints())
717 const auto& constraintInfo = constraintData.constraintInfo();
718 const auto& values = constraintData.values();
719 const auto dofIdx = constraintData.dofIndex();
721 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
723 if (constraintInfo.isConstraintEquation(eqIdx))
725 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
726 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
727 applyDirichletConstraint(dofIdx, values, eqIdx, pvIdx);
735 template<std::
size_t i,
class Res,
class GG,
class Sol>
736 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
738 if constexpr (Detail::hasPeriodicDofMap<GG>())
742 if (m.first < m.second)
745 res[m.first] += res[m.second];
748 res[m.second] = curSol[m.second] - curSol[m.first];
755 template<std::
size_t i,
class Res,
class GG,
class Sol>
756 void enforceProblemConstraints_(Dune::index_constant<i> domainI, Res& res,
const GG&
gridGeometry,
const Sol& curSol)
const
758 if constexpr (Detail::hasSubProblemGlobalConstraints<Problem<domainI>>())
760 for (
const auto& constraintData : this->
problem(domainI).constraints())
762 const auto& constraintInfo = constraintData.constraintInfo();
763 const auto& values = constraintData.values();
764 const auto dofIdx = constraintData.dofIndex();
765 for (
int eqIdx = 0; eqIdx < constraintInfo.size(); ++eqIdx)
767 if (constraintInfo.isConstraintEquation(eqIdx))
769 const auto pvIdx = constraintInfo.eqToPriVarIndex(eqIdx);
770 assert(0 <= pvIdx && pvIdx < constraintInfo.size());
771 res[dofIdx][eqIdx] = curSol[dofIdx][pvIdx] - values[pvIdx];
779 ProblemTuple problemTuple_;
782 GridGeometryTuple gridGeometryTuple_;
785 GridVariablesTuple gridVariablesTuple_;
788 std::shared_ptr<const TimeLoop> timeLoop_;
794 bool isStationaryProblem_;
797 std::shared_ptr<JacobianMatrix> jacobian_;
798 std::shared_ptr<ResidualType> residual_;
801 mutable bool warningIssued_;
804 bool enableMultithreading_ =
false;
Subdomain-specific views on multidomain assemblers.
Subdomain-specific view on a multidomain assembler. Allows retrieval of sub-domain specific objects w...
Definition: assemblerview.hh:31
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:103
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:380
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:361
void updateAfterGridAdaption()
Resizes jacobian and residual and recomputes colors.
Definition: multidomain/fvassembler.hh:351
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:125
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:122
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:413
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:133
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:277
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:124
CMType CouplingManager
Definition: multidomain/fvassembler.hh:128
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:270
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:371
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:405
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: multidomain/fvassembler.hh:435
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:116
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:390
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:319
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:395
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:119
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:409
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:428
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:421
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:400
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:453
MDTraits Traits
Definition: multidomain/fvassembler.hh:108
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager, std::shared_ptr< const TimeLoop > timeLoop, const SolutionVector &prevSol)
The constructor for instationary problems.
Definition: multidomain/fvassembler.hh:220
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:332
typename MDTraits::template SubDomain< id >::LocalResidual LocalResidual
Definition: multidomain/fvassembler.hh:113
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:448
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:385
typename MDTraits::ResidualVector ResidualType
Definition: multidomain/fvassembler.hh:126
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:192
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:249
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:441
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:110
ResidualType & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:417
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: multidomain/fvassembler.hh:304
The cell-centered scheme multidomain local assembler.
Definition: multidomain/subdomaincclocalassembler.hh:270
Definition: multidomain/fvassembler.hh:49
The default time loop for instationary simulations.
Definition: common/timeloop.hh:139
Defines all properties used in Dumux.
Manages the handling of time dependent problems.
Helper function to generate Jacobian pattern for multi domain models.
An enum class to define various differentiation methods available in order to compute the derivatives...
Some exceptions thrown in DuMux
dune-grid capabilities compatibility layer
constexpr bool isSerial()
Checking whether the backend is serial.
Definition: multithreading.hh:45
Helper function to generate Jacobian pattern for different discretization methods.
The available discretization methods in Dumux.
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain proble...
decltype(std::declval< P >().constraints()) SubProblemConstraintsDetector
helper struct detecting if sub-problem has a constraints() function
Definition: multidomain/fvassembler.hh:75
constexpr bool hasSubProblemGlobalConstraints()
Definition: multidomain/fvassembler.hh:78
bool allGridsSupportsMultithreadingImpl(const T &gridGeometries, std::index_sequence< I... >)
Definition: multistagemultidomainfvassembler.hh:52
bool allGridsSupportsMultithreading(const std::tuple< GG... > &gridGeometries)
Definition: multistagemultidomainfvassembler.hh:60
bool supportsMultithreading(const GridView &gridView)
Definition: gridcapabilities.hh:34
typename Detail::ConcatSeq< decltype(std::make_index_sequence< e >{}), e+1, decltype(std::make_index_sequence<(n > e) ?(n - e - 1) :0 >{})>::type makeIncompleteIntegerSequence
Definition: utility.hh:59
Provides a helper class for nonoverlapping decomposition.
Type traits to detect periodicity support.
Type trait that is specialized for coupling manager supporting multithreaded assembly.
Definition: multistagemultidomainfvassembler.hh:78
An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for...
Utilities for template meta programming.