13#ifndef DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
14#define DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
19#include <dune/common/indices.hh>
20#include <dune/common/exceptions.hh>
37template<
class MDTraits>
43 using Scalar =
typename MDTraits::Scalar;
44 using SolutionVector =
typename MDTraits::SolutionVector;
46 template<std::
size_t i>
using SubDomainTypeTag =
typename MDTraits::template SubDomain<i>::TypeTag;
52 template<std::
size_t i>
using GridGeometry =
typename MDTraits::template SubDomain<i>::GridGeometry;
53 template<std::
size_t i>
using FVElementGeometry =
typename GridGeometry<i>::LocalView;
54 template<std::
size_t i>
using SubControlVolumeFace =
typename GridGeometry<i>::SubControlVolumeFace;
55 template<std::
size_t i>
using SubControlVolume =
typename GridGeometry<i>::SubControlVolume;
56 template<std::
size_t i>
using GridView =
typename GridGeometry<i>::GridView;
57 template<std::
size_t i>
using Element =
typename GridView<i>::template Codim<0>::Entity;
59 template<std::
size_t i>
60 static constexpr auto domainIdx()
61 {
return typename MDTraits::template SubDomain<i>::Index{}; }
63 template<std::
size_t i>
64 static constexpr bool isCCTpfa()
67 using CouplingStencil = std::vector<std::size_t>;
79 void init(std::shared_ptr<Problem<domainIdx<0>()>> problem0,
80 std::shared_ptr<Problem<domainIdx<1>()>> problem1,
81 const SolutionVector&
curSol)
83 static_assert(isCCTpfa<0>() && isCCTpfa<1>(),
"Only cctpfa implemented!");
87 couplingMapper_.update(*
this);
105 template<std::
size_t i, std::
size_t j>
107 const Element<i>& element,
108 Dune::index_constant<j> domainJ)
const
110 static_assert(i != j,
"A domain cannot be coupled to itself!");
111 const auto eIdx = this->
problem(domainI).gridGeometry().elementMapper().index(element);
112 return couplingMapper_.couplingStencil(domainI, eIdx, domainJ);
122 template<std::
size_t i>
124 const SubControlVolumeFace<i>& scvf)
const
126 return couplingMapper_.isCoupled(domainI, scvf);
145 template<std::
size_t i>
147 const Element<i>& element,
148 const FVElementGeometry<i>& fvGeometry,
149 const ElementVolumeVariables<i>& elemVolVars,
150 const SubControlVolumeFace<i>& scvf,
151 int phaseIdx = 0)
const
154 static const bool enableGravity = getParamFromGroup<bool>(this->
problem(domainI).paramGroup(),
"Problem.EnableGravity");
155 constexpr auto otherDomainIdx = domainIdx<1-i>();
157 const auto& outsideElement = this->
problem(otherDomainIdx).gridGeometry().element(couplingMapper_.outsideElementIndex(domainI, scvf));
158 auto fvGeometryOutside =
localView(this->
problem(otherDomainIdx).gridGeometry());
159 fvGeometryOutside.bindElement(outsideElement);
161 const auto& flipScvf = fvGeometryOutside.scvf(couplingMapper_.flipScvfIndex(domainI, scvf));
162 const auto& outsideScv = fvGeometryOutside.scv(flipScvf.insideScvIdx());
163 const auto outsideVolVars =
volVars(otherDomainIdx, outsideElement, outsideScv);
165 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
166 const auto& insideVolVars = elemVolVars[insideScv];
168 using Extrusion =
typename GridGeometry<i>::Extrusion;
169 const auto ti =
computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideVolVars.permeability(), insideVolVars.extrusionFactor());
170 const auto tj =
computeTpfaTransmissibility(fvGeometryOutside, flipScvf, outsideScv, outsideVolVars.permeability(), outsideVolVars.extrusionFactor());
173 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
178 const auto rho = (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
179 const auto& g = this->
problem(domainI).spatialParams().gravity(scvf.ipGlobal());
182 const auto alpha_inside =
vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
183 const auto alpha_outside =
vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g)*outsideVolVars.extrusionFactor();
186 const auto pInside = insideVolVars.pressure(0);
187 const auto pOutside = outsideVolVars.pressure(0);
189 flux = tij*(pInside - pOutside) + rho*Extrusion::area(fvGeometry, scvf)*alpha_inside - rho*tij/tj*(alpha_inside - alpha_outside);
195 const auto pInside = insideVolVars.pressure(phaseIdx);
196 const auto pOutside = outsideVolVars.pressure(phaseIdx);
199 flux = tij*(pInside - pOutside);
203 static const Scalar upwindWeight = getParam<Scalar>(
"Flux.UpwindWeight");
204 auto upwindTerm = [phaseIdx](
const auto&
volVars){
return volVars.density(phaseIdx)*
volVars.mobility(phaseIdx); };
205 if (std::signbit(flux))
206 flux *= (upwindWeight*upwindTerm(outsideVolVars)
207 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
209 flux *= (upwindWeight*upwindTerm(insideVolVars)
210 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
213 flux /= Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor();
218 template<std::
size_t i>
219 VolumeVariables<i>
volVars(Dune::index_constant<i> domainI,
220 const Element<i>& element,
221 const SubControlVolume<i>& scv)
const
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:37
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:276
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:298
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:327
void updateSolution(const SolutionVector &curSol)
Updates the entire solution vector, e.g. before assembly or after grid adaption Overload might want t...
Definition: multidomain/couplingmanager.hh:208
Coupling manager for equal-dimension boundary coupling of darcy models.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:40
std::unordered_map< std::size_t, CouplingStencil > CouplingStencils
export stencil types
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:72
Scalar advectiveFluxCoupling(Dune::index_constant< i > domainI, const Element< i > &element, const FVElementGeometry< i > &fvGeometry, const ElementVolumeVariables< i > &elemVolVars, const SubControlVolumeFace< i > &scvf, int phaseIdx=0) const
Evaluate the boundary conditions for a coupled Neumann boundary segment.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:146
const CouplingStencil & couplingStencil(Dune::index_constant< i > domainI, const Element< i > &element, Dune::index_constant< j > domainJ) const
Methods to be accessed by the assembly.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:106
VolumeVariables< i > volVars(Dune::index_constant< i > domainI, const Element< i > &element, const SubControlVolume< i > &scv) const
Return the volume variables of domain i for a given element and scv.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:219
void init(std::shared_ptr< Problem< domainIdx< 0 >()> > problem0, std::shared_ptr< Problem< domainIdx< 1 >()> > problem1, const SolutionVector &curSol)
Methods to be accessed by main.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:79
bool isCoupled(Dune::index_constant< i > domainI, const SubControlVolumeFace< i > &scvf) const
If the boundary entity is on a coupling boundary.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:123
MDTraits MultiDomainTraits
export traits
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:70
the default mapper for conforming equal dimension boundary coupling between two domains (box or cc)
Definition: boundary/darcydarcy/couplingmapper.hh:37
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition: extrusion.hh:155
Defines all properties used in Dumux.
Element solution classes and factory functions.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:36
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:880
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:34
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
The interface of the coupling manager for multi domain problems.
constexpr CCTpfa cctpfa
Definition: method.hh:145
A helper to deduce a vector with the same size as numbers of equations.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...