3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
26#define DUMUX_MULTIDOMAIN_DARCYDARCY_BOUNDARY_COUPLINGMANAGER_HH
27
28#include <iostream>
29#include <vector>
30
31#include <dune/common/indices.hh>
32#include <dune/common/exceptions.hh>
33
35#include <dumux/common/math.hh>
42
43namespace Dumux {
44
49template<class MDTraits>
51: public CouplingManager<MDTraits>
52{
54
55 using Scalar = typename MDTraits::Scalar;
56 using SolutionVector = typename MDTraits::SolutionVector;
57
58 template<std::size_t i> using SubDomainTypeTag = typename MDTraits::template SubDomain<i>::TypeTag;
59 template<std::size_t i> using Problem = GetPropType<SubDomainTypeTag<i>, Properties::Problem>;
60 template<std::size_t i> using PrimaryVariables = GetPropType<SubDomainTypeTag<i>, Properties::PrimaryVariables>;
61 template<std::size_t i> using NumEqVector = Dumux::NumEqVector<PrimaryVariables<i>>;
62 template<std::size_t i> using ElementVolumeVariables = typename GetPropType<SubDomainTypeTag<i>, Properties::GridVolumeVariables>::LocalView;
63 template<std::size_t i> using VolumeVariables = typename GetPropType<SubDomainTypeTag<i>, Properties::GridVolumeVariables>::VolumeVariables;
64 template<std::size_t i> using GridGeometry = typename MDTraits::template SubDomain<i>::GridGeometry;
65 template<std::size_t i> using FVElementGeometry = typename GridGeometry<i>::LocalView;
66 template<std::size_t i> using SubControlVolumeFace = typename GridGeometry<i>::SubControlVolumeFace;
67 template<std::size_t i> using SubControlVolume = typename GridGeometry<i>::SubControlVolume;
68 template<std::size_t i> using GridView = typename GridGeometry<i>::GridView;
69 template<std::size_t i> using Element = typename GridView<i>::template Codim<0>::Entity;
70
71 template<std::size_t i>
72 static constexpr auto domainIdx()
73 { return typename MDTraits::template SubDomain<i>::Index{}; }
74
75 template<std::size_t i>
76 static constexpr bool isCCTpfa()
77 { return GridGeometry<i>::discMethod == DiscretizationMethods::cctpfa; }
78
79 using CouplingStencil = std::vector<std::size_t>;
80public:
82 using MultiDomainTraits = MDTraits;
84 using CouplingStencils = std::unordered_map<std::size_t, CouplingStencil>;
85
89 // \{
90
91 void init(std::shared_ptr<Problem<domainIdx<0>()>> problem0,
92 std::shared_ptr<Problem<domainIdx<1>()>> problem1,
93 const SolutionVector& curSol)
94 {
95 static_assert(isCCTpfa<0>() && isCCTpfa<1>(), "Only cctpfa implemented!");
96
97 this->setSubProblems(std::make_tuple(problem0, problem1));
98 this->updateSolution(curSol);
99 couplingMapper_.update(*this);
100 }
101
102 // \}
103
107 // \{
108
117 template<std::size_t i, std::size_t j>
118 const CouplingStencil& couplingStencil(Dune::index_constant<i> domainI,
119 const Element<i>& element,
120 Dune::index_constant<j> domainJ) const
121 {
122 static_assert(i != j, "A domain cannot be coupled to itself!");
123 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(element);
124 return couplingMapper_.couplingStencil(domainI, eIdx, domainJ);
125 }
126
127 // \}
128
134 template<std::size_t i>
135 bool isCoupled(Dune::index_constant<i> domainI,
136 const SubControlVolumeFace<i>& scvf) const
137 {
138 return couplingMapper_.isCoupled(domainI, scvf);
139 }
140
157 template<std::size_t i>
158 Scalar advectiveFluxCoupling(Dune::index_constant<i> domainI,
159 const Element<i>& element,
160 const FVElementGeometry<i>& fvGeometry,
161 const ElementVolumeVariables<i>& elemVolVars,
162 const SubControlVolumeFace<i>& scvf,
163 int phaseIdx = 0) const
164 {
165 Scalar flux = 0.0;
166 static const bool enableGravity = getParamFromGroup<bool>(this->problem(domainI).paramGroup(), "Problem.EnableGravity");
167 constexpr auto otherDomainIdx = domainIdx<1-i>();
168
169 const auto& outsideElement = this->problem(otherDomainIdx).gridGeometry().element(couplingMapper_.outsideElementIndex(domainI, scvf));
170 auto fvGeometryOutside = localView(this->problem(otherDomainIdx).gridGeometry());
171 fvGeometryOutside.bindElement(outsideElement);
172
173 const auto& flipScvf = fvGeometryOutside.scvf(couplingMapper_.flipScvfIndex(domainI, scvf));
174 const auto& outsideScv = fvGeometryOutside.scv(flipScvf.insideScvIdx());
175 const auto outsideVolVars = volVars(otherDomainIdx, outsideElement, outsideScv);
176
177 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
178 const auto& insideVolVars = elemVolVars[insideScv];
179
180 using Extrusion = typename GridGeometry<i>::Extrusion;
181 const auto ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, insideVolVars.permeability(), insideVolVars.extrusionFactor());
182 const auto tj = computeTpfaTransmissibility(fvGeometryOutside, flipScvf, outsideScv, outsideVolVars.permeability(), outsideVolVars.extrusionFactor());
183 Scalar tij = 0.0;
184 if (ti*tj > 0.0)
185 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
186
187 if (enableGravity)
188 {
189 // do averaging for the density over all neighboring elements
190 const auto rho = (insideVolVars.density(phaseIdx) + outsideVolVars.density(phaseIdx))*0.5;
191 const auto& g = this->problem(domainI).spatialParams().gravity(scvf.ipGlobal());
192
194 const auto alpha_inside = vtmv(scvf.unitOuterNormal(), insideVolVars.permeability(), g)*insideVolVars.extrusionFactor();
195 const auto alpha_outside = vtmv(scvf.unitOuterNormal(), outsideVolVars.permeability(), g)*outsideVolVars.extrusionFactor();
196
197 // Obtain inside and outside pressures
198 const auto pInside = insideVolVars.pressure(0);
199 const auto pOutside = outsideVolVars.pressure(0);
200
201 flux = tij*(pInside - pOutside) + rho*Extrusion::area(fvGeometry, scvf)*alpha_inside - rho*tij/tj*(alpha_inside - alpha_outside);
202
203 }
204 else
205 {
206 // Obtain inside and outside pressures
207 const auto pInside = insideVolVars.pressure(phaseIdx);
208 const auto pOutside = outsideVolVars.pressure(phaseIdx);
209
210 // return flux
211 flux = tij*(pInside - pOutside);
212 }
213
214 // upwind scheme
215 static const Scalar upwindWeight = getParam<Scalar>("Flux.UpwindWeight");
216 auto upwindTerm = [phaseIdx](const auto& volVars){ return volVars.density(phaseIdx)*volVars.mobility(phaseIdx); };
217 if (std::signbit(flux)) // if sign of flux is negative
218 flux *= (upwindWeight*upwindTerm(outsideVolVars)
219 + (1.0 - upwindWeight)*upwindTerm(insideVolVars));
220 else
221 flux *= (upwindWeight*upwindTerm(insideVolVars)
222 + (1.0 - upwindWeight)*upwindTerm(outsideVolVars));
223
224 // make this a area-specific flux
225 flux /= Extrusion::area(fvGeometry, scvf)*insideVolVars.extrusionFactor();
226 return flux;
227 }
228
230 template<std::size_t i>
231 VolumeVariables<i> volVars(Dune::index_constant<i> domainI,
232 const Element<i>& element,
233 const SubControlVolume<i>& scv) const
234 {
235 VolumeVariables<i> volVars;
236 const auto elemSol = elementSolution(element, this->curSol(domainI), this->problem(domainI).gridGeometry());
237 volVars.update(elemSol, this->problem(domainI), element, scv);
238 return volVars;
239 }
240
241private:
243};
244
245} // end namespace Dumux
246
247#endif
A helper to deduce a vector with the same size as numbers of equations.
Define some often used mathematical functions.
Element solution classes and factory functions.
The available discretization methods in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Tensor::field_type computeTpfaTransmissibility(const SubControlVolumeFace &scvf, const SubControlVolume &scv, const Tensor &T, typename 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:48
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:863
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:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr CCTpfa cctpfa
Definition: method.hh:134
A vector of primary variables.
Definition: common/properties.hh:49
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:55
The type for a global container for the volume variables.
Definition: common/properties.hh:107
Traits extracting the public Extrusion type from T Defaults to NoExtrusion if no such type is found.
Definition: extrusion.hh:240
Coupling manager for equal-dimension boundary coupling of darcy models.
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:52
std::unordered_map< std::size_t, CouplingStencil > CouplingStencils
export stencil types
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:84
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:158
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:118
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:231
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:91
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:135
MDTraits MultiDomainTraits
export traits
Definition: multidomain/boundary/darcydarcy/couplingmanager.hh:82
the default mapper for conforming equal dimension boundary coupling between two domains (box or cc)
Definition: boundary/darcydarcy/couplingmapper.hh:49
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
void setSubProblems(const std::tuple< std::shared_ptr< SubProblems >... > &problems)
set the pointers to the sub problems
Definition: multidomain/couplingmanager.hh:299
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
SubSolutionVector< i > & curSol(Dune::index_constant< i > domainIdx)
the solution vector of the subproblem
Definition: multidomain/couplingmanager.hh:350
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:231
Declares all properties used in Dumux.
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...
The interface of the coupling manager for multi domain problems.