3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggeredcouplingmanager.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 *****************************************************************************/
26#ifndef DUMUX_STAGGERED_COUPLING_MANAGER_HH
27#define DUMUX_STAGGERED_COUPLING_MANAGER_HH
28
33
34namespace Dumux {
35
41template<class MDTraits>
42class StaggeredCouplingManager: virtual public CouplingManager<MDTraits>
43{
45
46 template<std::size_t id> using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
47 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
48
49 using FVElementGeometry = typename GridGeometry<0>::LocalView;
50 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
51 using Element = typename GridView<0>::template Codim<0>::Entity;
52
53 using GridIndexType = typename IndexTraits< GridView<0> >::GridIndex;
54 using CouplingStencil = std::vector<GridIndexType>;
55
56public:
57
59
60 using Traits = MDTraits;
61
62 static constexpr auto cellCenterIdx = GridGeometry<0>::cellCenterIdx();
63 static constexpr auto faceIdx = GridGeometry<0>::faceIdx();
64
68 template<std::size_t i, std::size_t j, class LocalAssemblerI, class PriVarsJ>
69 void updateCouplingContext(Dune::index_constant<i> domainI,
70 const LocalAssemblerI& localAssemblerI,
71 Dune::index_constant<j> domainJ,
72 const std::size_t dofIdxGlobalJ,
73 const PriVarsJ& priVarsJ,
74 int pvIdxJ)
75 {
76 auto& curSol = this->curSol(domainJ);
77
78 // only proceed if the solution vector has actually been set in the base class
79 // (which is not the case if the staggered model is not coupled to another domain)
80 if (curSol.size() > 0)
81 curSol[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
82 }
83
94 const CouplingStencil& couplingStencil(Dune::index_constant<cellCenterIdx> domainI,
95 const Element& elementI,
96 Dune::index_constant<faceIdx> domainJ) const
97 {
98 const auto& connectivityMap = this->problem(domainI).gridGeometry().connectivityMap();
99 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(elementI);
100 return connectivityMap(domainI, domainJ, eIdx);
101 }
102
113 template<std::size_t i, std::size_t j>
114 const CouplingStencil couplingStencil(Dune::index_constant<i> domainI,
115 const SubControlVolumeFace& scvfI,
116 Dune::index_constant<j> domainJ) const
117 {
118 static_assert(i != j, "Domain i cannot be coupled to itself!");
119 static_assert(AlwaysFalse<Dune::index_constant<i>>::value,
120 "The coupling manager does not implement the couplingStencil() function" );
121
122 return CouplingStencil(); // supress compiler warning of function having no return statement
123 }
124
135 const CouplingStencil& couplingStencil(Dune::index_constant<faceIdx> domainI,
136 const SubControlVolumeFace& scvfI,
137 Dune::index_constant<cellCenterIdx> domainJ) const
138 {
139 const auto& connectivityMap = this->problem(domainI).gridGeometry().connectivityMap();
140 return connectivityMap(domainI, domainJ, scvfI.index());
141 }
142
149 template<class LocalAssemblerI, std::size_t j>
150 decltype(auto) evalCouplingResidual(Dune::index_constant<cellCenterIdx> domainI,
151 const LocalAssemblerI& localAssemblerI,
152 Dune::index_constant<j> domainJ,
153 std::size_t dofIdxGlobalJ) const
154 {
155 static_assert(domainI != domainJ, "Domain i cannot be coupled to itself!");
156 return localAssemblerI.evalLocalResidualForCellCenter();
157 }
158
176 template<class LocalAssemblerI, std::size_t j>
177 decltype(auto) evalCouplingResidual(Dune::index_constant<faceIdx> domainI,
178 const SubControlVolumeFace& scvfI,
179 const LocalAssemblerI& localAssemblerI,
180 Dune::index_constant<j> domainJ,
181 std::size_t dofIdxGlobalJ) const
182 {
183 static_assert(domainI != domainJ, "Domain i cannot be coupled to itself!");
184 return localAssemblerI.evalLocalResidualForFace(scvfI);
185 }
186
191 template<std::size_t i, typename std::enable_if_t<(GridGeometry<i>::discMethod != DiscretizationMethods::staggered), int> = 0>
192 decltype(auto) numericEpsilon(Dune::index_constant<i> id,
193 const std::string& paramGroup) const
194 {
195 return ParentType::numericEpsilon(id, paramGroup);
196 }
197
202 template<std::size_t i, typename std::enable_if_t<(GridGeometry<i>::discMethod == DiscretizationMethods::staggered), int> = 0>
203 decltype(auto) numericEpsilon(Dune::index_constant<i>,
204 const std::string& paramGroup) const
205 {
206 constexpr std::size_t numEqCellCenter = Traits::template SubDomain<cellCenterIdx>::PrimaryVariables::dimension;
207 constexpr std::size_t numEqFace = Traits::template SubDomain<faceIdx>::PrimaryVariables::dimension;
208 constexpr bool isCellCenter = GridGeometry<i>::isCellCenter();
209 constexpr std::size_t numEq = isCellCenter ? numEqCellCenter : numEqFace;
210 constexpr auto prefix = isCellCenter ? "CellCenter" : "Face";
211
212 try {
213 if(paramGroup.empty())
215 else
216 return NumericEpsilon<typename Traits::Scalar, numEq>(paramGroup + "." + prefix);
217 }
218 catch (Dune::RangeError& e)
219 {
220 DUNE_THROW(Dumux::ParameterException, "For the staggered model, you can to specify \n\n"
221 " CellCenter.Assembly.NumericDifference.PriVarMagnitude = mCC\n"
222 " Face.Assembly.NumericDifference.PriVarMagnitude = mFace\n"
223 " CellCenter.Assembly.NumericDifference.BaseEpsilon = eCC_0 ... eCC_numEqCellCenter-1\n"
224 " Face.Assembly.NumericDifference.BaseEpsilon = eFace_0 ... eFace_numEqFace-1\n\n"
225 "Wrong numer of values set for " << prefix << " (has " << numEq << " primary variable(s))\n\n" << e);
226 }
227
228 }
229
230};
231
232} //end namespace Dumux
233
234#endif
An adapter class for local assemblers using numeric differentiation.
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
decltype(auto) evalCouplingResidual(Dune::index_constant< faceIdx > domainI, const SubControlVolumeFace &scvfI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ) const
evaluates the face residual of a coupled face of domain i which depends on the variables at the degre...
Definition: staggeredcouplingmanager.hh:177
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...
Definition: multidomain/couplingmanager.hh:261
Definition: adapt.hh:29
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
Struture to define the index types used for grid and local indices.
Definition: indextraits.hh:38
Template which always yields a false value.
Definition: typetraits.hh:35
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
decltype(auto) curSol()
the solution vector of the coupled problem
Definition: multidomain/couplingmanager.hh:370
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:321
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
Definition: multidomain/couplingmanager.hh:287
Base coupling manager for the staggered discretization.
Definition: staggeredcouplingmanager.hh:43
const CouplingStencil & couplingStencil(Dune::index_constant< faceIdx > domainI, const SubControlVolumeFace &scvfI, Dune::index_constant< cellCenterIdx > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: staggeredcouplingmanager.hh:135
const CouplingStencil couplingStencil(Dune::index_constant< i > domainI, const SubControlVolumeFace &scvfI, Dune::index_constant< j > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: staggeredcouplingmanager.hh:114
decltype(auto) evalCouplingResidual(Dune::index_constant< cellCenterIdx > 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...
Definition: staggeredcouplingmanager.hh:150
const CouplingStencil & couplingStencil(Dune::index_constant< cellCenterIdx > domainI, const Element &elementI, Dune::index_constant< faceIdx > domainJ) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: staggeredcouplingmanager.hh:94
decltype(auto) numericEpsilon(Dune::index_constant< i > id, const std::string &paramGroup) const
return the numeric epsilon used for deflecting primary variables of coupled domain i.
Definition: staggeredcouplingmanager.hh:192
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, const std::size_t dofIdxGlobalJ, const PriVarsJ &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: staggeredcouplingmanager.hh:69
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.
Definition: staggeredcouplingmanager.hh:203
static constexpr auto cellCenterIdx
Definition: staggeredcouplingmanager.hh:62
MDTraits Traits
Definition: staggeredcouplingmanager.hh:60
static constexpr auto faceIdx
Definition: staggeredcouplingmanager.hh:63
The interface of the coupling manager for multi domain problems.