version 3.11-dev
Multidomain framework

Coupling of several regular DuMux problems. More...

Description

Multidomain Framework

Coupled models realized with the DuMux multidomain framework

DuMuX can couple problems posed on different domains. The domains can touch or overlap, model different physics, have different dimensions, different grids, or different discretization methods. The full system Jacobian is approximated by numeric differentiation (avoids hand-coding the Jacobian) which allows building monolithic solvers for complex nonlinear coupled problems. The framework enabling all these coupled simulations is called the MultiDomain framework.

Design Goals

Monolithic assembly and solvers

Structure of monolithic Jacobian matrix

The multidomain framework uses block-structured matrices to assemble a single system matrix for the linear system arising from the coupled problem (for example, in each Newton iteration of a nonlinear solver).

The matrix is structured such that each of the diagonal blocks correspond to one of the subproblems. The off-diagonal blocks correspond to the coupling derivatives that represent interactions between the subproblems. The block-structured matrix can directly be used with compatible linear solvers. The block structure is also useful to construct partitioned solvers.

The Subproblems

Any regular ("single domain") DuMux problem can be turned into a subproblem of a coupled simulation. The subproblems can model arbitrary physics, can have different dimensions, or be discretized with different methods. The only requirement is that the subproblems are compatible with the coupling manager.

In addition to a regular DuMux problem, the subproblems gets a pointer to the coupling manager which is used to transfer data between the subproblems and obtain data from other subproblems.

The Coupling Manager

The coupling manager is the central object in the multidomain framework, which is responsible for the coupling of the subproblems. Coupling managers implement the interface of the Dumux::CouplingManager class.

The coupling manager is responsible for the following tasks:


🛠 Edit above doc on GitLab

Modules

 Boundary coupling mode
 Couples problems of different or equal dimension that touch at the domain boundary. Examples are equal-dimension multi-physics problems like Darcy-Stokes coupling or PNM (pore network model)-Darcy coupling.
 
 Coupling for dual network approach for pore network models
 Coupling for dual network approach for pore network models.
 
 Embedded mixed-dimension coupling mode
 Couples problems of different dimensions where one or more lower-dimensional problems (lowdim) are embedded in a higher-dimensional domain (bulk). Examples are embedded one-dimensional networks for the simulation of blood tissue perfusion, or root-soil interaction, and embedded fracture models.
 
 Conforming mixed-dimension facet coupling mode
 Couples problems of different dimensions where one or more lower-dimensional problems (lowdim) live on the facets of the higher-dimensional domain (bulk). Examples are discrete facet conforming fracture models and problems with physics on a domain surface.
 

Classes

struct  Dumux::CouplingManagerSupportsMultithreadedAssembly< CM >
 Type trait that is specialized for coupling manager supporting multithreaded assembly. More...
 
class  Dumux::Experimental::MultiStageMultiDomainFVAssembler< MDTraits, CMType, diffMethod >
 A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa, mpfa, ...) with multiple domains. More...
 
class  Dumux::Experimental::SubDomainCCLocalAssemblerBase< id, TypeTag, Assembler, Implementation, dm >
 A base class for all multidomain local assemblers. More...
 
class  Dumux::Experimental::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DM >
 The cell-centered scheme multidomain local assembler. More...
 
class  Dumux::Experimental::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric >
 Cell-centered scheme multidomain local assembler using numeric differentiation. More...
 
class  Dumux::Experimental::SubDomainCVFELocalAssemblerBase< id, TypeTag, Assembler, Implementation, dm >
 A base class for all CVFE subdomain local assemblers. More...
 
class  Dumux::Experimental::SubDomainCVFELocalAssembler< id, TypeTag, Assembler, DM >
 The CVFE scheme multidomain local assembler. More...
 
class  Dumux::Experimental::SubDomainCVFELocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric >
 CVFE scheme multi domain local assembler using numeric differentiation. More...
 
class  Dumux::MultiDomainAssemblerSubDomainView< MDAssembler, domainId >
 Subdomain-specific view on a multidomain assembler. Allows retrieval of sub-domain specific objects w/o passing a domain id. More...
 
class  Dumux::CouplingManager< Traits >
 The interface of the coupling manager for multi domain problems. More...
 
class  Dumux::CVFEFreeFlowCouplingManager< Traits >
 The interface of the coupling manager for free flow systems. More...
 
class  Dumux::FCStaggeredFreeFlowCouplingManager< Traits >
 The interface of the coupling manager for free flow systems. More...
 
class  Dumux::MultiDomainFVAssembler< MDTraits, CMType, diffMethod, useImplicitAssembly >
 A linear system assembler (residual and Jacobian) for finite volume schemes (box, tpfa, mpfa, ...) with multiple domains. More...
 
class  Dumux::MultiDomainFVGridGeometry< MDTraits >
 A multidomain wrapper for multiple grid geometries. More...
 
class  Dumux::MultiDomainFVGridVariables< MDTraits >
 A multidomain wrapper for multiple grid variables. More...
 
class  Dumux::MultiDomainFVProblem< MDTraits >
 A multidomain wrapper for multiple problems. More...
 
class  Dumux::MultiDomainVtkOutputModule< MDTraits, Module >
 A multidomain wrapper for multiple vtk output modules. More...
 
class  Dumux::MultiBinaryCouplingManager< MDTraits, CouplingMap, CouplingMgrs >
 Coupling manager that combines an arbitrary number of binary coupling manager (coupling two domains each) More...
 
class  Dumux::MultiDomainNewtonConvergenceWriter< MDTraits >
 Writes the intermediate solutions for every Newton iteration. More...
 
class  Dumux::MultiDomainNewtonSolver< Assembler, LinearSolver, CouplingManager, Reassembler, Comm >
 Newton solver for coupled problems. More...
 
class  Dumux::SubDomainCCLocalAssemblerBase< id, TypeTag, Assembler, Implementation, implicit >
 A base class for all multidomain local assemblers. More...
 
class  Dumux::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DM, implicit >
 The cell-centered scheme multidomain local assembler. More...
 
class  Dumux::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, true >
 Cell-centered scheme multidomain local assembler using numeric differentiation and implicit time discretization. More...
 
class  Dumux::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, false >
 Cell-centered scheme multidomain local assembler using numeric differentiation and explicit time discretization. More...
 
class  Dumux::SubDomainCCLocalAssembler< id, TypeTag, Assembler, DiffMethod::analytic, true >
 Cell-centered scheme local assembler using analytic differentiation and implicit time discretization. More...
 
class  Dumux::SubDomainCVFELocalAssemblerBase< id, TypeTag, Assembler, Implementation, dm, implicit >
 A base class for all CVFE subdomain local assemblers. More...
 
class  Dumux::SubDomainCVFELocalAssembler< id, TypeTag, Assembler, DM, implicit >
 The CVFE scheme multidomain local assembler. More...
 
class  Dumux::SubDomainCVFELocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, true >
 CVFE scheme multi domain local assembler using numeric differentiation and implicit time discretization. More...
 
class  Dumux::SubDomainCVFELocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, false >
 CVFE scheme multi domain local assembler using numeric differentiation and explicit time discretization. More...
 
class  Dumux::SubDomainFaceCenteredLocalAssemblerBase< id, TypeTag, Assembler, Implementation, dm, implicit >
 A base class for all face-centered staggered local assemblers. More...
 
class  Dumux::SubDomainFaceCenteredLocalAssembler< id, TypeTag, Assembler, DM, implicit >
 The face-centered staggered scheme multidomain local assembler. More...
 
class  Dumux::SubDomainFaceCenteredLocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, true >
 Face-centered staggered scheme multi domain local assembler using numeric differentiation and implicit time discretization. More...
 
class  Dumux::SubDomainFaceCenteredLocalAssembler< id, TypeTag, Assembler, DiffMethod::numeric, false >
 Face-centered staggered scheme multi domain local assembler using numeric differentiation and explicit time discretization. More...
 

Typedefs

template<class Traits >
using Dumux::FreeFlowCouplingManager = typename Detail::FreeFlowCouplingManagerSelector< Traits >::type
 The interface of the coupling manager for free flow systems. More...
 
template<class DomainGridView , class TargetGridView , class DomainMapper , class TargetMapper >
using Dumux::MultiDomainGlue = IntersectionEntitySet< GridViewGeometricEntitySet< DomainGridView, 0, DomainMapper >, GridViewGeometricEntitySet< TargetGridView, 0, TargetMapper > >
 A convenience alias for the IntersectionEntitySet of two GridViewGeometricEntitySets. More...
 

Functions

template<bool isImplicit, class CouplingManager , class GridGeometryI , class GridGeometryJ , std::size_t i, std::size_t j, typename std::enable_if_t<((GridGeometryI::discMethod==DiscretizationMethods::cctpfa)||(GridGeometryI::discMethod==DiscretizationMethods::ccmpfa)), int > = 0>
Dune::MatrixIndexSet Dumux::getCouplingJacobianPattern (const CouplingManager &couplingManager, Dune::index_constant< i > domainI, const GridGeometryI &gridGeometryI, Dune::index_constant< j > domainJ, const GridGeometryJ &gridGeometryJ)
 Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for cell-centered schemes. More...
 
template<class DomainGG , class TargetGG >
MultiDomainGlue< typename DomainGG::GridView, typename TargetGG::GridView, typename DomainGG::ElementMapper, typename TargetGG::ElementMapper > Dumux::makeGlue (const DomainGG &domainGridGeometry, const TargetGG &targetGridGeometry)
 Creates the glue object containing the intersections between two grids obtained from given grid geometries. More...
 

Files

file  multistagemultidomainfvassembler.hh
 A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
 
file  experimental/assembly/subdomaincclocalassembler.hh
 A multidomain local assembler for Jacobian and residual contribution per element (cell-centered methods)
 
file  experimental/assembly/subdomaincvfelocalassembler.hh
 An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain problems.
 
file  assemblerview.hh
 Subdomain-specific views on multidomain assemblers.
 
file  couplingjacobianpattern.hh
 Helper function to generate Jacobian pattern for multi domain models.
 
file  multidomain/couplingmanager.hh
 The interface of the coupling manager for multi domain problems.
 
file  multidomain/freeflow/couplingmanager.hh
 Freeflow coupling managers (Navier-Stokes mass-momentum coupling)
 
file  couplingmanager_cvfe.hh
 Freeflow coupling managers (Navier-Stokes mass-momentum coupling)
 
file  couplingmanager_staggered.hh
 Freeflow coupling managers (Navier-Stokes mass-momentum coupling)
 
file  multidomain/freeflow/typetraits.hh
 Some useful type traits.
 
file  multidomain/fvassembler.hh
 A linear system assembler (residual and Jacobian) for finite volume schemes with multiple domains.
 
file  multidomain/fvgridgeometry.hh
 Multidomain wrapper for multiple grid geometries.
 
file  multidomain/fvgridvariables.hh
 Multidomain wrapper for multiple grid variables.
 
file  multidomain/fvproblem.hh
 Multidomain wrapper for multiple problems.
 
file  glue.hh
 A class glueing two grids of potentially different dimension geometrically. Intersections are computed using axis-aligned bounding box trees.
 
file  multidomain/io/vtkoutputmodule.hh
 Multidomain wrapper for multiple vtk output modules.
 
file  multibinarycouplingmanager.hh
 
file  multidomain/newtonconvergencewriter.hh
 This class provides the infrastructure to write the convergence behaviour of the Newton method for multidomain simulations into a VTK file.
 
file  multidomain/newtonsolver.hh
 
file  multidomain/subdomaincclocalassembler.hh
 A multidomain local assembler for Jacobian and residual contribution per element (cell-centered methods)
 
file  multidomain/subdomaincvfelocalassembler.hh
 An assembler for Jacobian and residual contribution per element (CVFE methods) for multidomain problems.
 
file  subdomainfclocalassembler.hh
 An assembler for Jacobian and residual contribution per element (face-centered staggered methods) for multidomain problems.
 
file  multidomain/traits.hh
 Traits for multidomain problems.
 

member functions concerning variable caching for element residual evaluations

template<std::size_t i, std::size_t j, class LocalAssemblerI >
void Dumux::CouplingManager< Traits >::updateCouplingContext (Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
 updates all data and variables that are necessary to evaluate the residual of the element of domain i this is called whenever one of the primary variables that the element residual depends on changes in domain j More...
 
template<std::size_t i, std::size_t j, class LocalAssemblerI >
decltype(auto) Dumux::CouplingManager< Traits >::evalCouplingResidual (Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
 evaluates the element residual of a coupled element of domain i which depends on the variables at the degree of freedom with index dofIdxGlobalJ of domain j More...
 

member functions concerning variable caching for element residual evaluations

template<std::size_t i, std::size_t j, class LocalAssemblerI >
void Dumux::CVFEFreeFlowCouplingManager< Traits >::updateCouplingContext (Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
 updates all data and variables that are necessary to evaluate the residual of the element of domain i this is called whenever one of the primary variables that the element residual depends on changes in domain j More...
 

Typedef Documentation

◆ FreeFlowCouplingManager

template<class Traits >
using Dumux::FreeFlowCouplingManager = typedef typename Detail::FreeFlowCouplingManagerSelector<Traits>::type

◆ MultiDomainGlue

template<class DomainGridView , class TargetGridView , class DomainMapper , class TargetMapper >
using Dumux::MultiDomainGlue = typedef IntersectionEntitySet<GridViewGeometricEntitySet<DomainGridView, 0, DomainMapper>, GridViewGeometricEntitySet<TargetGridView, 0, TargetMapper> >

Function Documentation

◆ evalCouplingResidual()

template<class Traits >
template<std::size_t i, std::size_t j, class LocalAssemblerI >
decltype(auto) Dumux::CouplingManager< Traits >::evalCouplingResidual ( Dune::index_constant< i >  domainI,
const LocalAssemblerI &  localAssemblerI,
Dune::index_constant< j >  domainJ,
std::size_t  dofIdxGlobalJ 
) const
inline
Parameters
domainIthe domain index of domain i
localAssemblerIthe local assembler assembling the element residual of an element of domain i
domainJthe domain index of domain j
dofIdxGlobalJthe index of the degree of freedom of domain j which has an influence on the element residual of domain i
Note
the element whose residual is to be evaluated can be retrieved from the local assembler as localAssemblerI.element() as well as all up-to-date variables and caches.
the default implementation evaluates the complete element residual if only parts (i.e. only certain scvs, or only certain terms of the residual) of the residual are coupled to dof with index dofIdxGlobalJ the function can be overloaded in the coupling manager
Returns
the element residual

◆ getCouplingJacobianPattern()

template<bool isImplicit, class CouplingManager , class GridGeometryI , class GridGeometryJ , std::size_t i, std::size_t j, typename std::enable_if_t<((GridGeometryI::discMethod==DiscretizationMethods::cctpfa)||(GridGeometryI::discMethod==DiscretizationMethods::ccmpfa)), int > = 0>
Dune::MatrixIndexSet Dumux::getCouplingJacobianPattern ( const CouplingManager couplingManager,
Dune::index_constant< i >  domainI,
const GridGeometryI &  gridGeometryI,
Dune::index_constant< j >  domainJ,
const GridGeometryJ &  gridGeometryJ 
)

Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for control-volume finite element schemes (CVFE)

Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for the staggered scheme (degrees of freedom on faces)

Helper function to generate coupling Jacobian pattern (off-diagonal blocks) for the staggered scheme (degrees of freedom on cell centers)

◆ makeGlue()

template<class DomainGG , class TargetGG >
MultiDomainGlue< typename DomainGG::GridView, typename TargetGG::GridView, typename DomainGG::ElementMapper, typename TargetGG::ElementMapper > Dumux::makeGlue ( const DomainGG &  domainGridGeometry,
const TargetGG &  targetGridGeometry 
)
Parameters
domainGridGeometryThe grid geometry of the domain
targetGridGeometryThe grid geometry of the target domain
Returns
The glue object containing the intersections

◆ updateCouplingContext() [1/2]

template<class Traits >
template<std::size_t i, std::size_t j, class LocalAssemblerI >
void Dumux::CouplingManager< Traits >::updateCouplingContext ( Dune::index_constant< i >  domainI,
const LocalAssemblerI &  localAssemblerI,
Dune::index_constant< j >  domainJ,
std::size_t  dofIdxGlobalJ,
const PrimaryVariables< j > &  priVarsJ,
int  pvIdxJ 
)
inline
Parameters
domainIthe domain index of domain i
localAssemblerIthe local assembler assembling the element residual of an element of domain i
domainJthe domain index of domain j
dofIdxGlobalJthe index of the degree of freedom of domain j whose solution changed
priVarsJthe new solution at the degree of freedom of domain j with index dofIdxGlobalJ
pvIdxJthe index of the primary variable of domain j which has been updated
Note
this concerns all data that is used in the evaluation of the element residual and depends on the primary variables at the degree of freedom location with index dofIdxGlobalJ
the element whose residual is to be evaluated can be retrieved from the local assembler as localAssemblerI.element()
per default, we update the solution vector, if the element residual of domain i depends on more than the primary variables of domain j update the other dependent data here by overloading this function

◆ updateCouplingContext() [2/2]

template<class Traits >
template<std::size_t i, std::size_t j, class LocalAssemblerI >
void Dumux::CVFEFreeFlowCouplingManager< Traits >::updateCouplingContext ( Dune::index_constant< i >  domainI,
const LocalAssemblerI &  localAssemblerI,
Dune::index_constant< j >  domainJ,
std::size_t  dofIdxGlobalJ,
const PrimaryVariables< j > &  priVarsJ,
int  pvIdxJ 
)
inline
Parameters
domainIthe domain index of domain i
localAssemblerIthe local assembler assembling the element residual of an element of domain i
domainJthe domain index of domain j
dofIdxGlobalJthe index of the degree of freedom of domain j whose solution changed
priVarsJthe new solution at the degree of freedom of domain j with index dofIdxGlobalJ
pvIdxJthe index of the primary variable of domain j which has been updated
Note
this concerns all data that is used in the evaluation of the element residual and depends on the primary variables at the degree of freedom location with index dofIdxGlobalJ
the element whose residual is to be evaluated can be retrieved from the local assembler as localAssemblerI.element()
per default, we update the solution vector, if the element residual of domain i depends on more than the primary variables of domain j update the other dependent data here by overloading this function