3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Classes | Public Types | Public Member Functions | Static Public Attributes | List of all members
Dumux::FCStaggeredFreeFlowCouplingManager< Traits > Class Template Reference

The interface of the coupling manager for free flow systems. More...

#include <dumux/multidomain/freeflow/couplingmanager_staggered.hh>

Inheritance diagram for Dumux::FCStaggeredFreeFlowCouplingManager< Traits >:

Description

template<class Traits>
class Dumux::FCStaggeredFreeFlowCouplingManager< Traits >

The interface of the coupling manager for free flow systems.

Note
coupling manager the face-centered staggered discretization scheme

Public Types

using SolutionVectorStorage = typename ParentType::SolutionVectorStorage
 

Public Member Functions

void init (std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol)
 Methods to be accessed by main. More...
 
void init (std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, const SolutionVector &curSol, const SolutionVector &prevSol)
 use as regular coupling manager in a transient setting More...
 
void init (std::shared_ptr< Problem< freeFlowMomentumIndex > > momentumProblem, std::shared_ptr< Problem< freeFlowMassIndex > > massProblem, GridVariablesTuple &&gridVariables, typename ParentType::SolutionVectorStorage &curSol)
 use as binary coupling manager in multi model context More...
 
template<std::size_t j, class LocalAssemblerI >
decltype(auto) evalCouplingResidual (Dune::index_constant< freeFlowMomentumIndex > domainI, const LocalAssemblerI &localAssemblerI, const SubControlVolume< freeFlowMomentumIndex > &scvI, 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...
 
template<std::size_t i, std::size_t j, class LocalAssemblerI >
decltype(auto) 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 the coupling stencils
Scalar pressure (const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
 Returns the pressure at a given frontal sub control volume face. More...
 
Scalar cellPressure (const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
 Returns the pressure at the center of a sub control volume corresponding to a given sub control volume face. This is used for setting a Dirichlet pressure for the mass model when a fixed pressure for the momentum balance is set at another boundary. Since the the pressure at the given scvf is solution-dependent and thus unknown a priori, we just use the value of the interior cell here. More...
 
Scalar density (const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
 Returns the density at a given sub control volume face. More...
 
auto insideAndOutsideDensity (const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf, const bool considerPreviousTimeStep=false) const
 
Scalar density (const Element< freeFlowMomentumIndex > &element, const SubControlVolume< freeFlowMomentumIndex > &scv, const bool considerPreviousTimeStep=false) const
 Returns the density at a given sub control volume. More...
 
Scalar effectiveViscosity (const Element< freeFlowMomentumIndex > &element, const FVElementGeometry< freeFlowMomentumIndex > &fvGeometry, const SubControlVolumeFace< freeFlowMomentumIndex > &scvf) const
 Returns the pressure at a given sub control volume face. More...
 
VelocityVector faceVelocity (const Element< freeFlowMassIndex > &element, const SubControlVolumeFace< freeFlowMassIndex > &scvf) const
 Returns the velocity at a given sub control volume face. More...
 
template<std::size_t j>
const CouplingStencilTypecouplingStencil (Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< j > domainJ) const
 The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depends on. More...
 
const CouplingStencilTypecouplingStencil (Dune::index_constant< freeFlowMassIndex > domainI, const Element< freeFlowMassIndex > &elementI, Dune::index_constant< freeFlowMomentumIndex > domainJ) const
 returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the element residual of the given element of domain i More...
 
const CouplingStencilTypecouplingStencil (Dune::index_constant< freeFlowMomentumIndex > domainI, const Element< freeFlowMomentumIndex > &elementI, const SubControlVolume< freeFlowMomentumIndex > &scvI, Dune::index_constant< freeFlowMassIndex > domainJ) const
 returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the residual of the given sub-control volume of domain i More...
 
member functions concerning the coupling stencils
template<std::size_t i, std::size_t j>
const CouplingStencilType< i, j > & couplingStencil (Dune::index_constant< i > domainI, const Element< i > &elementI, Dune::index_constant< j > domainJ) const
 returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the element residual of the given element of domain i More...
 
template<std::size_t id, class JacobianPattern >
void extendJacobianPattern (Dune::index_constant< id > domainI, JacobianPattern &pattern) const
 extend the jacobian pattern of the diagonal block of domain i by those entries that are not already in the uncoupled pattern More...
 

Static Public Attributes

static constexpr auto freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index()
 
static constexpr auto freeFlowMassIndex = typename Traits::template SubDomain<1>::Index()
 
static constexpr auto pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx
 

member functions concerning variable caching for element residual evaluations

template<std::size_t i, std::size_t j, class LocalAssemblerI >
void 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...
 
void computeColorsForAssembly ()
 Compute colors for multithreaded assembly. More...
 
template<std::size_t i, class AssembleElementFunc >
void assembleMultithreaded (Dune::index_constant< i > domainI, AssembleElementFunc &&assembleElement) const
 Execute assembly kernel in parallel. More...
 

member functions concerning variable caching for element residual evaluations

template<std::size_t i, class Assembler >
void bindCouplingContext (Dune::index_constant< i > domainI, const Element< i > &elementI, const Assembler &assembler)
 prepares all data and variables that are necessary to evaluate the residual of the element of domain i More...
 
template<std::size_t i, class LocalAssemblerI , class UpdatableElementVolVars , class UpdatableFluxVarCache >
void updateCoupledVariables (Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, UpdatableElementVolVars &elemVolVars, UpdatableFluxVarCache &elemFluxVarsCache)
 update variables of domain i that depend on variables in domain j after the coupling context has been updated More...
 
void updateSolution (const SolutionVector &curSol)
 Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want to overload function if the solution vector is stored outside this class to make sure updates don't happen more than once. More...
 
template<std::size_t i, class LocalAssemblerI , class JacobianMatrixDiagBlock , class GridVariables >
void evalAdditionalDomainDerivatives (Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const typename LocalAssemblerI::LocalResidual::ElementResidualVector &origResiduals, JacobianMatrixDiagBlock &A, GridVariables &gridVariables)
 evaluate additional derivatives of the element residual of a domain with respect to dofs in the same domain that are not in the regular stencil (see CouplingManager::extendJacobianPattern) More...
 
template<std::size_t i>
decltype(auto) numericEpsilon (Dune::index_constant< i >, const std::string &paramGroup) const
 return the numeric epsilon used for deflecting primary variables of coupled domain i More...
 
template<typename... SubProblems>
void setSubProblems (const std::tuple< std::shared_ptr< SubProblems >... > &problems)
 set the pointers to the sub problems More...
 
template<class SubProblem , std::size_t i>
void setSubProblem (std::shared_ptr< SubProblem > problem, Dune::index_constant< i > domainIdx)
 set a pointer to one of the sub problems More...
 
template<std::size_t i>
const Problem< i > & problem (Dune::index_constant< i > domainIdx) const
 Return a reference to the sub problem. More...
 
void attachSolution (SolutionVectorStorage &curSol)
 Attach a solution vector stored outside of this class. More...
 
template<std::size_t i>
SubSolutionVector< i > & curSol (Dune::index_constant< i > domainIdx)
 the solution vector of the subproblem More...
 
template<std::size_t i>
const SubSolutionVector< i > & curSol (Dune::index_constant< i > domainIdx) const
 the solution vector of the subproblem More...
 

Member Typedef Documentation

◆ SolutionVectorStorage

template<class Traits >
using Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::SolutionVectorStorage = typename ParentType::SolutionVectorStorage

Member Function Documentation

◆ assembleMultithreaded()

template<class Traits >
template<std::size_t i, class AssembleElementFunc >
void Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::assembleMultithreaded ( Dune::index_constant< i >  domainI,
AssembleElementFunc &&  assembleElement 
) const
inline

Execute assembly kernel in parallel.

Parameters
domainIthe domain index of domain i
assembleElementkernel function to execute for one element

◆ attachSolution()

template<class Traits >
void Dumux::CouplingManager< Traits >::attachSolution ( SolutionVectorStorage curSol)
inlineprotectedinherited

Attach a solution vector stored outside of this class.

Note
The caller has to make sure that curSol stays alive for the lifetime of the coupling manager. Otherwise we have a dangling reference here. Use with care.

◆ bindCouplingContext()

template<class Traits >
template<std::size_t i, class Assembler >
void Dumux::CouplingManager< Traits >::bindCouplingContext ( Dune::index_constant< i >  domainI,
const Element< i > &  elementI,
const Assembler &  assembler 
)
inlineinherited

prepares all data and variables that are necessary to evaluate the residual of the element of domain i

Parameters
domainIthe domain index of domain i
elementIthe element whose residual we are assemling next
assemblerthe multidomain assembler for access to all data necessary for the assembly of all domains
Note
this concerns all data that is used in the evaluation of the element residual and depends on one of the degrees of freedom returned by CouplingManager::couplingStencil
every coupled element residual depends at least on the solution of another domain, that why we always store a copy of the solution vector in the coupling manager, hence, in case the element residual only depends on primary variables of the other domain this function does nothing
overload this function in case the element residual depends on more than the primary variables of domain j

◆ cellPressure()

template<class Traits >
Scalar Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::cellPressure ( const Element< freeFlowMassIndex > &  element,
const SubControlVolumeFace< freeFlowMassIndex > &  scvf 
) const
inline

Returns the pressure at the center of a sub control volume corresponding to a given sub control volume face. This is used for setting a Dirichlet pressure for the mass model when a fixed pressure for the momentum balance is set at another boundary. Since the the pressure at the given scvf is solution-dependent and thus unknown a priori, we just use the value of the interior cell here.

◆ computeColorsForAssembly()

template<class Traits >
void Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::computeColorsForAssembly ( )
inline

Compute colors for multithreaded assembly.

Parameters
domainIthe domain index of domain i
assembleElementkernel function to execute for one element

◆ couplingStencil() [1/4]

template<class Traits >
const CouplingStencilType & Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::couplingStencil ( Dune::index_constant< freeFlowMassIndex domainI,
const Element< freeFlowMassIndex > &  elementI,
Dune::index_constant< freeFlowMomentumIndex domainJ 
) const
inline

returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the element residual of the given element of domain i

Parameters
domainIthe domain index of domain i
elementIthe coupled element of domain í
domainJthe domain index of domain j
Note
The element residual definition depends on the discretization scheme of domain i box: a container of the residuals of all sub control volumes cc : the residual of the (sub) control volume fem: the residual of the element
This function has to be implemented by all coupling managers for all combinations of i and j

◆ couplingStencil() [2/4]

template<class Traits >
const CouplingStencilType & Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::couplingStencil ( Dune::index_constant< freeFlowMomentumIndex domainI,
const Element< freeFlowMomentumIndex > &  elementI,
const SubControlVolume< freeFlowMomentumIndex > &  scvI,
Dune::index_constant< freeFlowMassIndex domainJ 
) const
inline

returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the residual of the given sub-control volume of domain i

Parameters
domainIthe domain index of domain i
elementIthe coupled element of domain í
scvIthe sub-control volume of domain i
domainJthe domain index of domain j

◆ couplingStencil() [3/4]

template<class Traits >
template<std::size_t j>
const CouplingStencilType & Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::couplingStencil ( Dune::index_constant< freeFlowMomentumIndex domainI,
const Element< freeFlowMomentumIndex > &  elementI,
const SubControlVolume< freeFlowMomentumIndex > &  scvI,
Dune::index_constant< j >  domainJ 
) const
inline

The coupling stencil of domain I, i.e. which domain J DOFs the given domain I element's residual depends on.

◆ couplingStencil() [4/4]

template<class Traits >
template<std::size_t i, std::size_t j>
const CouplingStencilType< i, j > & Dumux::CouplingManager< Traits >::couplingStencil ( Dune::index_constant< i >  domainI,
const Element< i > &  elementI,
Dune::index_constant< j >  domainJ 
) const
inlineinherited

returns an iterable container of all indices of degrees of freedom of domain j that couple with / influence the element residual of the given element of domain i

Parameters
domainIthe domain index of domain i
elementIthe coupled element of domain í
domainJthe domain index of domain j
Note
The element residual definition depends on the discretization scheme of domain i box: a container of the residuals of all sub control volumes cc : the residual of the (sub) control volume fem: the residual of the element
This function has to be implemented by all coupling managers for all combinations of i and j

◆ curSol() [1/2]

template<class Traits >
template<std::size_t i>
SubSolutionVector< i > & Dumux::CouplingManager< Traits >::curSol ( Dune::index_constant< i >  domainIdx)
inlineprotectedinherited

the solution vector of the subproblem

Parameters
domainIdxThe domain index
Note
in case of numeric differentiation the solution vector always carries the deflected solution

◆ curSol() [2/2]

template<class Traits >
template<std::size_t i>
const SubSolutionVector< i > & Dumux::CouplingManager< Traits >::curSol ( Dune::index_constant< i >  domainIdx) const
inlineprotectedinherited

the solution vector of the subproblem

Parameters
domainIdxThe domain index
Note
in case of numeric differentiation the solution vector always carries the deflected solution

◆ density() [1/2]

template<class Traits >
Scalar Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::density ( const Element< freeFlowMomentumIndex > &  element,
const FVElementGeometry< freeFlowMomentumIndex > &  fvGeometry,
const SubControlVolumeFace< freeFlowMomentumIndex > &  scvf,
const bool  considerPreviousTimeStep = false 
) const
inline

Returns the density at a given sub control volume face.

◆ density() [2/2]

template<class Traits >
Scalar Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::density ( const Element< freeFlowMomentumIndex > &  element,
const SubControlVolume< freeFlowMomentumIndex > &  scv,
const bool  considerPreviousTimeStep = false 
) const
inline

Returns the density at a given sub control volume.

◆ effectiveViscosity()

template<class Traits >
Scalar Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::effectiveViscosity ( const Element< freeFlowMomentumIndex > &  element,
const FVElementGeometry< freeFlowMomentumIndex > &  fvGeometry,
const SubControlVolumeFace< freeFlowMomentumIndex > &  scvf 
) const
inline

Returns the pressure at a given sub control volume face.

◆ evalAdditionalDomainDerivatives()

template<class Traits >
template<std::size_t i, class LocalAssemblerI , class JacobianMatrixDiagBlock , class GridVariables >
void Dumux::CouplingManager< Traits >::evalAdditionalDomainDerivatives ( Dune::index_constant< i >  domainI,
const LocalAssemblerI &  localAssemblerI,
const typename LocalAssemblerI::LocalResidual::ElementResidualVector &  origResiduals,
JacobianMatrixDiagBlock &  A,
GridVariables &  gridVariables 
)
inlineinherited

evaluate additional derivatives of the element residual of a domain with respect to dofs in the same domain that are not in the regular stencil (see CouplingManager::extendJacobianPattern)

Note
Such additional dependencies can arise from the coupling, e.g. if a coupling source term depends on a non-local average of a quantity of the same domain

◆ 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

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

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

◆ extendJacobianPattern()

template<class Traits >
template<std::size_t id, class JacobianPattern >
void Dumux::CouplingManager< Traits >::extendJacobianPattern ( Dune::index_constant< id >  domainI,
JacobianPattern &  pattern 
) const
inlineinherited

extend the jacobian pattern of the diagonal block of domain i by those entries that are not already in the uncoupled pattern

Note
per default we do not add such additional dependencies
Such additional dependencies can arise from the coupling, e.g. if a coupling source term depends on a non-local average of a quantity of the same domain
Warning List:
if you overload this also implement evalAdditionalDomainDerivatives

◆ faceVelocity()

template<class Traits >
VelocityVector Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::faceVelocity ( const Element< freeFlowMassIndex > &  element,
const SubControlVolumeFace< freeFlowMassIndex > &  scvf 
) const
inline

Returns the velocity at a given sub control volume face.

◆ init() [1/3]

template<class Traits >
void Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::init ( std::shared_ptr< Problem< freeFlowMomentumIndex > >  momentumProblem,
std::shared_ptr< Problem< freeFlowMassIndex > >  massProblem,
GridVariablesTuple &&  gridVariables,
const SolutionVector &  curSol 
)
inline

Methods to be accessed by main.

use as regular coupling manager

◆ init() [2/3]

template<class Traits >
void Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::init ( std::shared_ptr< Problem< freeFlowMomentumIndex > >  momentumProblem,
std::shared_ptr< Problem< freeFlowMassIndex > >  massProblem,
GridVariablesTuple &&  gridVariables,
const SolutionVector &  curSol,
const SolutionVector &  prevSol 
)
inline

use as regular coupling manager in a transient setting

◆ init() [3/3]

template<class Traits >
void Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::init ( std::shared_ptr< Problem< freeFlowMomentumIndex > >  momentumProblem,
std::shared_ptr< Problem< freeFlowMassIndex > >  massProblem,
GridVariablesTuple &&  gridVariables,
typename ParentType::SolutionVectorStorage curSol 
)
inline

use as binary coupling manager in multi model context

◆ insideAndOutsideDensity()

template<class Traits >
auto Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::insideAndOutsideDensity ( const Element< freeFlowMomentumIndex > &  element,
const FVElementGeometry< freeFlowMomentumIndex > &  fvGeometry,
const SubControlVolumeFace< freeFlowMomentumIndex > &  scvf,
const bool  considerPreviousTimeStep = false 
) const
inline

◆ numericEpsilon()

template<class Traits >
template<std::size_t i>
decltype(auto) Dumux::CouplingManager< Traits >::numericEpsilon ( Dune::index_constant< i >  ,
const std::string &  paramGroup 
) const
inlineinherited

return the numeric epsilon used for deflecting primary variables of coupled domain i

◆ pressure()

template<class Traits >
Scalar Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::pressure ( const Element< freeFlowMomentumIndex > &  element,
const FVElementGeometry< freeFlowMomentumIndex > &  fvGeometry,
const SubControlVolumeFace< freeFlowMomentumIndex > &  scvf 
) const
inline

Returns the pressure at a given frontal sub control volume face.

◆ problem()

template<class Traits >
template<std::size_t i>
const Problem< i > & Dumux::CouplingManager< Traits >::problem ( Dune::index_constant< i >  domainIdx) const
inlineinherited

Return a reference to the sub problem.

Parameters
domainIdxThe domain index We avoid exception handling here because the performance of this function is critical

◆ setSubProblem()

template<class Traits >
template<class SubProblem , std::size_t i>
void Dumux::CouplingManager< Traits >::setSubProblem ( std::shared_ptr< SubProblem >  problem,
Dune::index_constant< i >  domainIdx 
)
inlineinherited

set a pointer to one of the sub problems

Parameters
problema pointer to the sub problem
domainIdxthe domain index of the sub problem

◆ setSubProblems()

template<class Traits >
template<typename... SubProblems>
void Dumux::CouplingManager< Traits >::setSubProblems ( const std::tuple< std::shared_ptr< SubProblems >... > &  problems)
inlineinherited

set the pointers to the sub problems

Parameters
problemsA tuple of shared pointers to the sub problems

◆ updateCoupledVariables()

template<class Traits >
template<std::size_t i, class LocalAssemblerI , class UpdatableElementVolVars , class UpdatableFluxVarCache >
void Dumux::CouplingManager< Traits >::updateCoupledVariables ( Dune::index_constant< i >  domainI,
const LocalAssemblerI &  localAssemblerI,
UpdatableElementVolVars &  elemVolVars,
UpdatableFluxVarCache &  elemFluxVarsCache 
)
inlineinherited

update variables of domain i that depend on variables in domain j after the coupling context has been updated

Parameters
domainIthe index of domain i
localAssemblerIthe local assembler assembling the element residual of an element of domain i
elemVolVarsthe element volume variables (all volume variables in the element local stencil) to be updated
elemFluxVarsCachethe element flux variable cache (all flux variables in the element local stencil) to be updated
Note
Such variables do not necessarily exist and then this function does nothing (default)
some examples from geomechanics: the porosity of (physical) domain i (porous medium flow) depends on the displacement vector of physical domain j (mechanics) from domaindecomposition: the transmissibilities for fluxes of domain i to domain j depend on the permeability in domain j (which might depend in turn on the primary variables of domain i)

◆ updateCouplingContext()

template<class Traits >
template<std::size_t i, std::size_t j, class LocalAssemblerI >
void Dumux::FCStaggeredFreeFlowCouplingManager< 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

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

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

◆ updateSolution()

template<class Traits >
void Dumux::CouplingManager< Traits >::updateSolution ( const SolutionVector curSol)
inlineinherited

Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want to overload function if the solution vector is stored outside this class to make sure updates don't happen more than once.

Member Data Documentation

◆ freeFlowMassIndex

template<class Traits >
constexpr auto Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::freeFlowMassIndex = typename Traits::template SubDomain<1>::Index()
staticconstexpr

◆ freeFlowMomentumIndex

template<class Traits >
constexpr auto Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::freeFlowMomentumIndex = typename Traits::template SubDomain<0>::Index()
staticconstexpr

◆ pressureIdx

template<class Traits >
constexpr auto Dumux::FCStaggeredFreeFlowCouplingManager< Traits >::pressureIdx = VolumeVariables<freeFlowMassIndex>::Indices::pressureIdx
staticconstexpr

The documentation for this class was generated from the following file: