26#ifndef DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
27#define DUMUX_MULTIDOMAIN_FV_ASSEMBLER_HH
32#include <dune/common/hybridutilities.hh>
33#include <dune/istl/matrixindexset.hh>
62template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
65 template<std::
size_t id>
66 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
71 using Scalar =
typename MDTraits::Scalar;
74 template<std::
size_t id>
77 template<std::
size_t id>
78 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
80 template<std::
size_t id>
81 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
83 template<std::
size_t id>
84 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
96 {
return useImplicitAssembly; }
100 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
101 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
102 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
107 template<
class DiscretizationMethod, std::
size_t id>
108 struct SubDomainAssemblerType;
110 template<std::
size_t id>
111 struct SubDomainAssemblerType<DiscretizationMethods::CCTpfa, id>
116 template<std::
size_t id>
117 struct SubDomainAssemblerType<DiscretizationMethods::CCMpfa, id>
119 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
122 template<std::
size_t id>
123 struct SubDomainAssemblerType<DiscretizationMethods::Box, id>
125 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
128 template<std::
size_t id>
129 struct SubDomainAssemblerType<DiscretizationMethods::Staggered, id>
131 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
134 template<std::
size_t id>
135 struct SubDomainAssemblerType<DiscretizationMethods::FCStaggered, id>
137 using type = SubDomainFaceCenteredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
140 template<std::
size_t id>
156 , problemTuple_(std::move(
problem))
160 , isStationaryProblem_(true)
161 , warningIssued_(false)
163 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
164 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
176 std::shared_ptr<const TimeLoop> timeLoop,
179 , problemTuple_(std::move(
problem))
182 , timeLoop_(timeLoop)
184 , isStationaryProblem_(false)
185 , warningIssued_(false)
187 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
196 checkAssemblerState_();
200 using namespace Dune::Hybrid;
201 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
203 auto& jacRow = (*jacobian_)[domainId];
204 auto& subRes = (*residual_)[domainId];
205 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
207 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
208 enforcePeriodicConstraints_(domainId, jacRow, subRes, *
gridGeometry, curSol[domainId]);
224 checkAssemblerState_();
229 using namespace Dune::Hybrid;
230 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
232 auto& subRes = r[domainId];
233 this->assembleResidual_(domainId, subRes, curSol);
245 Scalar resultSquared = 0.0;
248 using namespace Dune::Hybrid;
249 forEach(integralRange(Dune::Hybrid::size(
residual)), [&](
const auto domainId)
251 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
264 vectorHelper.makeNonOverlappingConsistent(
residual[domainId]);
267 else if (!warningIssued_)
270 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
271 <<
"overlapping entities multiple times. Please use the norm\n"
272 <<
"function provided by a linear solver instead." << std::endl;
274 warningIssued_ =
true;
281 localNormSquared =
gridView.comm().sum(localNormSquared);
284 resultSquared += localNormSquared;
288 return sqrt(resultSquared);
297 std::shared_ptr<SolutionVector> r)
313 jacobian_ = std::make_shared<JacobianMatrix>();
314 residual_ = std::make_shared<SolutionVector>();
326 using namespace Dune::Hybrid;
327 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
329 forEach(jac[i], [&](
auto& jacBlock)
331 using BlockType = std::decay_t<
decltype(jacBlock)>;
332 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
333 jacBlock.setBuildMode(BlockType::BuildMode::random);
334 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
335 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
345 using namespace Dune::Hybrid;
346 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
348 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
350 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
351 pattern.exportIdx(jac[domainI][domainJ]);
361 using namespace Dune::Hybrid;
362 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
363 { res[domainId].resize(this->
numDofs(domainId)); });
371 using namespace Dune::Hybrid;
372 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
373 { this->
gridVariables(domainId).update(curSol[domainId]); });
381 using namespace Dune::Hybrid;
382 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
383 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
387 template<std::
size_t i>
388 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
389 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
392 template<std::
size_t i>
393 const auto&
problem(Dune::index_constant<i> domainId)
const
394 {
return *std::get<domainId>(problemTuple_); }
397 template<std::
size_t i>
399 {
return *std::get<domainId>(gridGeometryTuple_); }
402 template<std::
size_t i>
403 const auto&
gridView(Dune::index_constant<i> domainId)
const
407 template<std::
size_t i>
409 {
return *std::get<domainId>(gridVariablesTuple_); }
412 template<std::
size_t i>
414 {
return *std::get<domainId>(gridVariablesTuple_); }
422 {
return *jacobian_; }
426 {
return *residual_; }
430 {
return *prevSol_; }
437 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
450 {
return isStationaryProblem_; }
455 template<std::
size_t i>
457 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
465 void resetResidual_()
469 residual_ = std::make_shared<SolutionVector>();
477 void resetJacobian_()
481 jacobian_ = std::make_shared<JacobianMatrix>();
490 void checkAssemblerState_()
const
492 if (!isStationaryProblem_ && !prevSol_)
493 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
495 if (isStationaryProblem_ && prevSol_)
496 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
497 <<
" Did you forget to set the timeLoop to make this problem instationary?");
500 template<std::
size_t i,
class JacRow,
class SubRes>
501 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
504 assemble_(domainId, [&](
const auto& element)
506 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
507 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
511 template<std::
size_t i,
class SubRes>
512 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
515 assemble_(domainId, [&](
const auto& element)
517 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
518 subDomainAssembler.assembleResidual(subRes);
528 template<std::
size_t i,
class AssembleElementFunc>
529 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
532 for (
const auto& element : elements(
gridView(domainId)))
533 assembleElement(element);
537 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
538 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
539 Dune::index_constant<j> domainJ)
const
542 auto pattern = getJacobianPattern<isImplicit()>(gg);
548 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
549 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
550 Dune::index_constant<j> domainJ)
const
558 template<std::
size_t i,
class JacRow,
class Sol,
class GG>
559 void enforcePeriodicConstraints_(Dune::index_constant<i> domainI, JacRow& jacRow, Sol& res,
const GG&
gridGeometry,
const Sol& curSol)
565 if (m.first < m.second)
567 auto& jac = jacRow[domainI];
570 res[m.first] += res[m.second];
573 res[m.second] = curSol[m.second] - curSol[m.first];
575 const auto end = jac[m.second].end();
576 for (
auto it = jac[m.second].begin(); it != end; ++it)
577 jac[m.first][it.index()] += (*it);
580 for (
auto it = jac[m.second].begin(); it != end; ++it)
581 (*it) = it.index() == m.second ? 1.0 : it.index() == m.first ? -1.0 : 0.0;
583 using namespace Dune::Hybrid;
586 auto& jacCoupling = jacRow[couplingDomainId];
588 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
589 jacCoupling[m.first][it.index()] += (*it);
591 for (
auto it = jacCoupling[m.second].begin(); it != jacCoupling[m.second].end(); ++it)
600 ProblemTuple problemTuple_;
603 GridGeometryTuple gridGeometryTuple_;
606 GridVariablesTuple gridVariablesTuple_;
609 std::shared_ptr<const TimeLoop> timeLoop_;
615 bool isStationaryProblem_;
618 std::shared_ptr<JacobianMatrix> jacobian_;
619 std::shared_ptr<SolutionVector> residual_;
An enum class to define various differentiation methods available in order to compute the derivatives...
Helper function to generate Jacobian pattern for different discretization methods.
Utilities for template meta programming.
The available discretization methods in Dumux.
Provides a helper class for nonoverlapping decomposition.
Helper function to generate Jacobian pattern for multi domain models.
An assembler for Jacobian and residual contribution per element (box methods) for multidomain problem...
A multidomain local assembler for Jacobian and residual contribution per element (cell-centered metho...
An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for...
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:42
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:71
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
typename BlockTypeHelper< SolutionVector, Dune::IsNumber< SolutionVector >::value >::type BlockType
Definition: nonlinear/newtonsolver.hh:198
constexpr Box box
Definition: method.hh:139
constexpr FCStaggered fcstaggered
Definition: method.hh:142
Definition: common/properties.hh:74
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
Definition: parallelhelpers.hh:430
A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa,...
Definition: multidomain/fvassembler.hh:64
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:388
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:369
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:87
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:88
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:84
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:421
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:95
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:220
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:86
CMType CouplingManager
Definition: multidomain/fvassembler.hh:90
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:213
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:379
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:238
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:413
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:443
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:78
void setLinearSystem(std::shared_ptr< JacobianMatrix > A, std::shared_ptr< SolutionVector > r)
Tells the assembler which jacobian and residual to use. This also resizes the containers to the requi...
Definition: multidomain/fvassembler.hh:296
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:398
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:359
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:311
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:403
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:81
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:417
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:436
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:429
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:408
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:461
MDTraits Traits
Definition: multidomain/fvassembler.hh:69
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:172
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:324
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:75
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:456
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:393
MultiDomainFVAssembler(ProblemTuple problem, GridGeometryTuple gridGeometry, GridVariablesTuple gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:151
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:194
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:449
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:71
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:343
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:425
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Declares all properties used in Dumux.
Manages the handling of time dependent problems.