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>
61template<
class MDTraits,
class CMType, DiffMethod diffMethod,
bool useImplicitAssembly = true>
64 template<std::
size_t id>
65 using SubDomainTypeTag =
typename MDTraits::template SubDomain<id>::TypeTag;
70 using Scalar =
typename MDTraits::Scalar;
73 template<std::
size_t id>
76 template<std::
size_t id>
77 using GridVariables =
typename MDTraits::template SubDomain<id>::GridVariables;
79 template<std::
size_t id>
80 using GridGeometry =
typename MDTraits::template SubDomain<id>::GridGeometry;
82 template<std::
size_t id>
83 using Problem =
typename MDTraits::template SubDomain<id>::Problem;
95 {
return useImplicitAssembly; }
99 using ProblemTuple =
typename MDTraits::template TupleOfSharedPtrConst<Problem>;
100 using GridGeometryTuple =
typename MDTraits::template TupleOfSharedPtrConst<GridGeometry>;
101 using GridVariablesTuple =
typename MDTraits::template TupleOfSharedPtr<GridVariables>;
106 template<DiscretizationMethod discMethod, std::
size_t id>
107 struct SubDomainAssemblerType;
109 template<std::
size_t id>
115 template<std::
size_t id>
118 using type = SubDomainCCLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
121 template<std::
size_t id>
124 using type = SubDomainBoxLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
127 template<std::
size_t id>
130 using type = SubDomainStaggeredLocalAssembler<id, SubDomainTypeTag<id>, ThisType, diffMethod,
isImplicit()>;
133 template<std::
size_t id>
134 using SubDomainAssembler =
typename SubDomainAssemblerType<GridGeometry<id>::discMethod,
id>::type;
153 , isStationaryProblem_(true)
154 , warningIssued_(false)
156 static_assert(
isImplicit(),
"Explicit assembler for stationary problem doesn't make sense!");
157 std::cout <<
"Instantiated assembler for a stationary problem." << std::endl;
169 std::shared_ptr<const TimeLoop> timeLoop,
175 , timeLoop_(timeLoop)
177 , isStationaryProblem_(false)
178 , warningIssued_(false)
180 std::cout <<
"Instantiated assembler for an instationary problem." << std::endl;
189 checkAssemblerState_();
193 using namespace Dune::Hybrid;
194 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainId)
196 auto& jacRow = (*jacobian_)[domainId];
197 auto& subRes = (*residual_)[domainId];
198 this->assembleJacobianAndResidual_(domainId, jacRow, subRes, curSol);
214 checkAssemblerState_();
219 using namespace Dune::Hybrid;
220 forEach(integralRange(Dune::Hybrid::size(r)), [&](
const auto domainId)
222 auto& subRes = r[domainId];
223 this->assembleResidual_(domainId, subRes, curSol);
235 Scalar resultSquared = 0.0;
238 using namespace Dune::Hybrid;
239 forEach(integralRange(Dune::Hybrid::size(
residual)), [&](
const auto domainId)
241 const auto gridGeometry = std::get<domainId>(gridGeometryTuple_);
254 vectorHelper.makeNonOverlappingConsistent(
residual[domainId]);
257 else if (!warningIssued_)
260 std::cout <<
"\nWarning: norm calculation adds entries corresponding to\n"
261 <<
"overlapping entities multiple times. Please use the norm\n"
262 <<
"function provided by a linear solver instead." << std::endl;
264 warningIssued_ =
true;
271 localNormSquared =
gridView.comm().sum(localNormSquared);
274 resultSquared += localNormSquared;
278 return sqrt(resultSquared);
287 std::shared_ptr<SolutionVector> r)
303 jacobian_ = std::make_shared<JacobianMatrix>();
304 residual_ = std::make_shared<SolutionVector>();
316 using namespace Dune::Hybrid;
317 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto i)
319 forEach(jac[i], [&](
auto& jacBlock)
321 using BlockType = std::decay_t<
decltype(jacBlock)>;
322 if (jacBlock.buildMode() == BlockType::BuildMode::unknown)
323 jacBlock.setBuildMode(BlockType::BuildMode::random);
324 else if (jacBlock.buildMode() != BlockType::BuildMode::random)
325 DUNE_THROW(Dune::NotImplemented,
"Only BCRS matrices with random build mode are supported at the moment");
335 using namespace Dune::Hybrid;
336 forEach(std::make_index_sequence<JacobianMatrix::N()>(), [&](
const auto domainI)
338 forEach(integralRange(Dune::Hybrid::size(jac[domainI])), [&](
const auto domainJ)
340 const auto pattern = this->getJacobianPattern_(domainI, domainJ);
341 pattern.exportIdx(jac[domainI][domainJ]);
351 using namespace Dune::Hybrid;
352 forEach(integralRange(Dune::Hybrid::size(res)), [&](
const auto domainId)
353 { res[domainId].resize(this->
numDofs(domainId)); });
361 using namespace Dune::Hybrid;
362 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
363 { this->
gridVariables(domainId).update(curSol[domainId]); });
371 using namespace Dune::Hybrid;
372 forEach(integralRange(Dune::Hybrid::size(gridVariablesTuple_)), [&](
const auto domainId)
373 { this->
gridVariables(domainId).resetTimeStep(curSol[domainId]); });
377 template<std::
size_t i>
378 std::size_t
numDofs(Dune::index_constant<i> domainId)
const
379 {
return std::get<domainId>(gridGeometryTuple_)->numDofs(); }
382 template<std::
size_t i>
383 const auto&
problem(Dune::index_constant<i> domainId)
const
384 {
return *std::get<domainId>(problemTuple_); }
387 template<std::
size_t i>
389 {
return *std::get<domainId>(gridGeometryTuple_); }
392 template<std::
size_t i>
393 const auto&
gridView(Dune::index_constant<i> domainId)
const
397 template<std::
size_t i>
399 {
return *std::get<domainId>(gridVariablesTuple_); }
402 template<std::
size_t i>
404 {
return *std::get<domainId>(gridVariablesTuple_); }
412 {
return *jacobian_; }
416 {
return *residual_; }
420 {
return *prevSol_; }
427 { timeLoop_ = timeLoop; isStationaryProblem_ = !(
static_cast<bool>(timeLoop)); }
440 {
return isStationaryProblem_; }
445 template<std::
size_t i>
447 {
return LocalResidual<i>(std::get<domainId>(problemTuple_).get(), timeLoop_.get()); }
455 void resetResidual_()
459 residual_ = std::make_shared<SolutionVector>();
467 void resetJacobian_()
471 jacobian_ = std::make_shared<JacobianMatrix>();
480 void checkAssemblerState_()
const
482 if (!isStationaryProblem_ && !prevSol_)
483 DUNE_THROW(Dune::InvalidStateException,
"Assembling instationary problem but previous solution was not set!");
485 if (isStationaryProblem_ && prevSol_)
486 DUNE_THROW(Dune::InvalidStateException,
"Assembling stationary problem but a previous solution was set."
487 <<
" Did you forget to set the timeLoop to make this problem instationary?");
490 template<std::
size_t i,
class JacRow,
class SubRes>
491 void assembleJacobianAndResidual_(Dune::index_constant<i> domainId, JacRow& jacRow, SubRes& subRes,
494 assemble_(domainId, [&](
const auto& element)
496 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
497 subDomainAssembler.assembleJacobianAndResidual(jacRow, subRes, gridVariablesTuple_);
501 template<std::
size_t i,
class SubRes>
502 void assembleResidual_(Dune::index_constant<i> domainId, SubRes& subRes,
505 assemble_(domainId, [&](
const auto& element)
507 SubDomainAssembler<i> subDomainAssembler(*
this, element, curSol, *
couplingManager_);
508 subDomainAssembler.assembleResidual(subRes);
518 template<std::
size_t i,
class AssembleElementFunc>
519 void assemble_(Dune::index_constant<i> domainId, AssembleElementFunc&& assembleElement)
const
522 for (
const auto& element : elements(
gridView(domainId)))
523 assembleElement(element);
527 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i==j),
int> = 0>
528 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
529 Dune::index_constant<j> domainJ)
const
532 auto pattern = getJacobianPattern<isImplicit()>(gg);
538 template<std::
size_t i, std::
size_t j,
typename std::enable_if_t<(i!=j),
int> = 0>
539 Dune::MatrixIndexSet getJacobianPattern_(Dune::index_constant<i> domainI,
540 Dune::index_constant<j> domainJ)
const
548 ProblemTuple problemTuple_;
551 GridGeometryTuple gridGeometryTuple_;
554 GridVariablesTuple gridVariablesTuple_;
557 std::shared_ptr<const TimeLoop> timeLoop_;
563 bool isStationaryProblem_;
566 std::shared_ptr<JacobianMatrix> jacobian_;
567 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.
A multidomain assembler for Jacobian and residual contribution per element (staggered method)
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...
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
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:63
std::size_t numDofs(Dune::index_constant< i > domainId) const
the number of dof locations of domain i
Definition: multidomain/fvassembler.hh:378
void updateGridVariables(const SolutionVector &curSol)
Updates the grid variables with the given solution.
Definition: multidomain/fvassembler.hh:359
typename MDTraits::SolutionVector SolutionVector
Definition: multidomain/fvassembler.hh:86
SolutionVector ResidualType
Definition: multidomain/fvassembler.hh:87
typename MDTraits::template SubDomain< id >::Problem Problem
Definition: multidomain/fvassembler.hh:83
JacobianMatrix & jacobian()
the full Jacobian matrix
Definition: multidomain/fvassembler.hh:411
static constexpr bool isImplicit()
Returns true if the assembler considers implicit assembly.
Definition: multidomain/fvassembler.hh:94
void assembleResidual(ResidualType &r, const SolutionVector &curSol)
assemble a residual r
Definition: multidomain/fvassembler.hh:210
typename MDTraits::JacobianMatrix JacobianMatrix
Definition: multidomain/fvassembler.hh:85
CMType CouplingManager
Definition: multidomain/fvassembler.hh:89
void assembleResidual(const SolutionVector &curSol)
compute the residuals using the internal residual
Definition: multidomain/fvassembler.hh:203
void resetTimeStep(const SolutionVector &curSol)
Resets the grid variables to the last time step.
Definition: multidomain/fvassembler.hh:369
Scalar residualNorm(const SolutionVector &curSol)
compute the residual and return it's vector norm
Definition: multidomain/fvassembler.hh:228
const GridVariables< i > & gridVariables(Dune::index_constant< i > domainId) const
the grid variables of domain i
Definition: multidomain/fvassembler.hh:403
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:433
typename MDTraits::template SubDomain< id >::GridVariables GridVariables
Definition: multidomain/fvassembler.hh:77
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:286
const auto & gridGeometry(Dune::index_constant< i > domainId) const
the finite volume grid geometry of domain i
Definition: multidomain/fvassembler.hh:388
void setResidualSize(SolutionVector &res) const
Resizes the residual.
Definition: multidomain/fvassembler.hh:349
void setLinearSystem()
The version without arguments uses the default constructor to create the jacobian and residual object...
Definition: multidomain/fvassembler.hh:301
const auto & gridView(Dune::index_constant< i > domainId) const
the grid view of domain i
Definition: multidomain/fvassembler.hh:393
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:165
typename MDTraits::template SubDomain< id >::GridGeometry GridGeometry
Definition: multidomain/fvassembler.hh:80
const CouplingManager & couplingManager() const
the coupling manager
Definition: multidomain/fvassembler.hh:407
void setTimeManager(std::shared_ptr< const TimeLoop > timeLoop)
Set time loop for instationary problems.
Definition: multidomain/fvassembler.hh:426
const SolutionVector & prevSol() const
the solution of the previous time step
Definition: multidomain/fvassembler.hh:419
GridVariables< i > & gridVariables(Dune::index_constant< i > domainId)
the grid variables of domain i
Definition: multidomain/fvassembler.hh:398
std::shared_ptr< CouplingManager > couplingManager_
the coupling manager coupling the sub domains
Definition: multidomain/fvassembler.hh:451
MDTraits Traits
Definition: multidomain/fvassembler.hh:68
void setJacobianBuildMode(JacobianMatrix &jac) const
Sets the jacobian build mode.
Definition: multidomain/fvassembler.hh:314
GetPropType< SubDomainTypeTag< id >, Properties::LocalResidual > LocalResidual
TODO get rid of this GetPropType.
Definition: multidomain/fvassembler.hh:74
MultiDomainFVAssembler(ProblemTuple &&problem, GridGeometryTuple &&gridGeometry, GridVariablesTuple &&gridVariables, std::shared_ptr< CouplingManager > couplingManager)
The constructor for stationary problems.
Definition: multidomain/fvassembler.hh:144
LocalResidual< i > localResidual(Dune::index_constant< i > domainId) const
Create a local residual object (used by the local assembler)
Definition: multidomain/fvassembler.hh:446
const auto & problem(Dune::index_constant< i > domainId) const
the problem of domain i
Definition: multidomain/fvassembler.hh:383
void assembleJacobianAndResidual(const SolutionVector &curSol)
Assembles the global Jacobian of the residual and the residual for the current solution.
Definition: multidomain/fvassembler.hh:187
bool isStationaryProblem() const
Whether we are assembling a stationary or instationary problem.
Definition: multidomain/fvassembler.hh:439
typename MDTraits::Scalar Scalar
Definition: multidomain/fvassembler.hh:70
void setJacobianPattern(JacobianMatrix &jac) const
Sets the jacobian sparsity pattern.
Definition: multidomain/fvassembler.hh:333
SolutionVector & residual()
the full residual vector
Definition: multidomain/fvassembler.hh:415
The cell-centered scheme multidomain local assembler.
Definition: subdomaincclocalassembler.hh:282
Declares all properties used in Dumux.
Manages the handling of time dependent problems.