version 3.10-dev
multidomain/boundary/darcydarcy/couplingmanager.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
14#define DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
15
16#include <iostream>
17#include <vector>
18
19#include <dune/common/indices.hh>
20#include <dune/common/exceptions.hh>
21
23#include <dumux/common/math.hh>
30
31namespace Dumux {
32
37template<class MDTraits>
39: public CouplingManager<MDTraits>
40{
42
43 using Scalar = typename MDTraits::Scalar;
44 using SolutionVector = typename MDTraits::SolutionVector;
45
46 template<std::size_t i> using SubDomainTypeTag = typename MDTraits::template SubDomain<i>::TypeTag;
47 template<std::size_t i> using Problem = GetPropType<SubDomainTypeTag<i>, Properties::Problem>;
48 template<std::size_t i> using PrimaryVariables = GetPropType<SubDomainTypeTag<i>, Properties::PrimaryVariables>;
49 template<std::size_t i> using NumEqVector = Dumux::NumEqVector<PrimaryVariables<i>>;
50 template<std::size_t i> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<i>, Properties::GridVolumeVariables>::LocalView;
51 template<std::size_t i> using VolumeVariables = typename GetPropType<SubDomainTypeTag<i>, Properties::GridVolumeVariables>::VolumeVariables;
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;
58
59 template<std::size_t i>
60 static constexpr auto domainIdx()
61 { return typename MDTraits::template SubDomain<i>::Index{}; }
62
63 template<std::size_t i>
64 static constexpr bool isCCTpfa()
65 { return GridGeometry<i>::discMethod == DiscretizationMethods::cctpfa; }
66
67 using CouplingStencil = std::vector<std::size_t>;
68public:
70 using MultiDomainTraits = MDTraits;
72 using CouplingStencils = std::unordered_map<std::size_t, CouplingStencil>;
73
77 // \{
78
79 void init(std::shared_ptr<Problem<domainIdx<0>()>> problem0,
80 std::shared_ptr<Problem<domainIdx<1>()>> problem1,
81 const SolutionVector& curSol)
82 {
83 static_assert(isCCTpfa<0>() && isCCTpfa<1>(), "Only cctpfa implemented!");
84
85 this->setSubProblems(std::make_tuple(problem0, problem1));
86 this->updateSolution(curSol);
87 couplingMapper_.update(*this);
88 }
89
90 // \}
91
95 // \{
96
105 template<std::size_t i, std::size_t j>
106 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
107 const Element<i>& element,
108 Dune::index_constant<j> domainJ) const
109 {
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);
113 }
114
115 // \}
116
122 template<std::size_t i>
123 bool isCoupled(Dune::index_constant<i> domainI,
124 const SubControlVolumeFace<i>& scvf) const
125 {
126 return couplingMapper_.isCoupled(domainI, scvf);
127 }
128
145 template<std::size_t i>
146 Scalar advectiveFluxCoupling(Dune::index_constant<i> domainI,
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
152 {
153 Scalar flux = 0.0;
154 static const bool enableGravity = getParamFromGroup<bool>(this->problem(domainI).paramGroup(), "Problem.EnableGravity");
155 constexpr auto otherDomainIdx = domainIdx<1-i>();
156
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);
160
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);
164
165 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
166 const auto& insideVolVars = elemVolVars[insideScv];
167
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());
171 Scalar tij = 0.0;
172 if (ti*tj > 0.0)
173 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
174
175 if (enableGravity)
176 {
177 // do averaging for the density over all neighboring elements
178 const auto rho = (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
179 const auto& g = this->problem(domainI).spatialParams().gravity(scvf.ipGlobal());
180
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();
184
185 // Obtain inside and outside pressures
186 const auto pInside = insideVolVars.pressure(0);
187 const auto pOutside = outsideVolVars.pressure(0);
188
189 flux = tij*(pInside - pOutside) + rho*Extrusion::area(fvGeometry, scvf)*alpha_inside - rho*tij/tj*(alpha_inside - alpha_outside);
190
191 }
192 else
193 {
194 // Obtain inside and outside pressures
195 const auto pInside = insideVolVars.pressure(phaseIdx);
196 const auto pOutside = outsideVolVars.pressure(phaseIdx);
197
198 // return flux
199 flux = tij*(pInside - pOutside);
200 }
201
202 // upwind scheme
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)) // if sign of flux is negative
206 flux *= (upwindWeight*upwindTerm(outsideVolVars)
207 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
208 else
209 flux *= (upwindWeight*upwindTerm(insideVolVars)
210 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
211
212 // make this a area-specific flux
213 flux /= Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor();
214 return flux;
215 }
216
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
222 {
223 VolumeVariables<i> volVars;
224 const auto elemSol = elementSolution(element, this->curSol(domainI), this->problem(domainI).gridGeometry());
225 volVars.update(elemSol, this->problem(domainI), element, scv);
226 return volVars;
227 }
228
229private:
231};
232
233} // end namespace Dumux
234
235#endif
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
Definition: adapt.hh:17
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...