GCC Code Coverage Report


Directory: ../../../builds/dumux-repositories/
File: dumux/dumux/multidomain/staggeredcouplingmanager.hh
Date: 2025-04-12 19:19:20
Exec Total Coverage
Lines: 27 29 93.1%
Functions: 65 65 100.0%
Branches: 40 48 83.3%

Line Branch Exec Source
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5 // SPDX-License-Identifier: GPL-3.0-or-later
6 //
7 /*!
8 * \file
9 * \ingroup MultiDomain
10 * \ingroup StaggeredDiscretization
11 * \brief The interface of the coupling manager for multi domain problems
12 */
13
14 #ifndef DUMUX_STAGGERED_COUPLING_MANAGER_HH
15 #define DUMUX_STAGGERED_COUPLING_MANAGER_HH
16
17 #include <dumux/multidomain/couplingmanager.hh>
18 #include <dumux/assembly/numericepsilon.hh>
19 #include <dumux/common/indextraits.hh>
20 #include <dumux/common/typetraits/typetraits.hh>
21 #include <dumux/discretization/method.hh>
22
23 namespace Dumux {
24
25 /*!
26 * \ingroup MultiDomain
27 * \ingroup StaggeredDiscretization
28 * \brief Base coupling manager for the staggered discretization.
29 */
30 template<class MDTraits>
31 85 class StaggeredCouplingManager: virtual public CouplingManager<MDTraits>
32 {
33 using ParentType = CouplingManager<MDTraits>;
34
35 template<std::size_t id> using GridGeometry = typename MDTraits::template SubDomain<id>::GridGeometry;
36 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
37
38 using FVElementGeometry = typename GridGeometry<0>::LocalView;
39 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
40 using Element = typename GridView<0>::template Codim<0>::Entity;
41
42 using GridIndexType = typename IndexTraits< GridView<0> >::GridIndex;
43 using CouplingStencil = std::vector<GridIndexType>;
44
45 public:
46
47 using ParentType::evalCouplingResidual;
48
49 using Traits = MDTraits;
50
51 static constexpr auto cellCenterIdx = GridGeometry<0>::cellCenterIdx();
52 static constexpr auto faceIdx = GridGeometry<0>::faceIdx();
53
54 /*!
55 * \copydoc Dumux::CouplingManager::updateCouplingContext()
56 */
57 template<std::size_t i, std::size_t j, class LocalAssemblerI, class PriVarsJ>
58 293469106 void updateCouplingContext(Dune::index_constant<i> domainI,
59 const LocalAssemblerI& localAssemblerI,
60 Dune::index_constant<j> domainJ,
61 const std::size_t dofIdxGlobalJ,
62 const PriVarsJ& priVarsJ,
63 int pvIdxJ)
64 {
65
16/16
✓ Branch 0 taken 6983200 times.
✓ Branch 1 taken 34757732 times.
✓ Branch 2 taken 6983200 times.
✓ Branch 3 taken 68787264 times.
✓ Branch 4 taken 1652160 times.
✓ Branch 5 taken 3790128 times.
✓ Branch 6 taken 1652160 times.
✓ Branch 7 taken 7334256 times.
✓ Branch 8 taken 2893408 times.
✓ Branch 9 taken 28380556 times.
✓ Branch 10 taken 2893408 times.
✓ Branch 11 taken 14341578 times.
✓ Branch 12 taken 13016320 times.
✓ Branch 13 taken 57179744 times.
✓ Branch 14 taken 13016320 times.
✓ Branch 15 taken 29807672 times.
293469106 auto& curSol = this->curSol(domainJ);
66
67 // only proceed if the solution vector has actually been set in the base class
68 // (which is not the case if the staggered model is not coupled to another domain)
69
16/16
✓ Branch 0 taken 6983200 times.
✓ Branch 1 taken 34757732 times.
✓ Branch 2 taken 6983200 times.
✓ Branch 3 taken 68787264 times.
✓ Branch 4 taken 1652160 times.
✓ Branch 5 taken 3790128 times.
✓ Branch 6 taken 1652160 times.
✓ Branch 7 taken 7334256 times.
✓ Branch 8 taken 2893408 times.
✓ Branch 9 taken 28380556 times.
✓ Branch 10 taken 2893408 times.
✓ Branch 11 taken 14341578 times.
✓ Branch 12 taken 13016320 times.
✓ Branch 13 taken 57179744 times.
✓ Branch 14 taken 13016320 times.
✓ Branch 15 taken 29807672 times.
293469106 if (curSol.size() > 0)
70 49090176 curSol[dofIdxGlobalJ][pvIdxJ] = priVarsJ[pvIdxJ];
71 }
72
73 /*!
74 * \brief returns an iterable container of all indices of degrees of freedom of domain j
75 * that couple with / influence the element residual of the given element of domain i
76 *
77 * \param domainI the domain index of domain i
78 * \param elementI the coupled element of domain í
79 * \param domainJ the domain index of domain j
80
81 * \note this is a specialization for getting the indices of the coupled staggered face dofs
82 */
83 89592 const CouplingStencil& couplingStencil(Dune::index_constant<cellCenterIdx> domainI,
84 const Element& elementI,
85 Dune::index_constant<faceIdx> domainJ) const
86 {
87
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 89592 times.
89592 const auto& connectivityMap = this->problem(domainI).gridGeometry().connectivityMap();
88 89592 const auto eIdx = this->problem(domainI).gridGeometry().elementMapper().index(elementI);
89 89592 return connectivityMap(domainI, domainJ, eIdx);
90 }
91
92 /*!
93 * \brief returns an iterable container of all indices of degrees of freedom of domain j
94 * that couple with / influence the residual of the given subcontrolvolume face of domain i
95 *
96 * \param domainI the domain index of domain i
97 * \param scvfI the coupled subcontrolvolume face of domain í
98 * \param domainJ the domain index of domain j
99
100 * \note This function has to be implemented by all coupling managers for all combinations of i and j
101 */
102 template<std::size_t i, std::size_t j>
103 const CouplingStencil couplingStencil(Dune::index_constant<i> domainI,
104 const SubControlVolumeFace& scvfI,
105 Dune::index_constant<j> domainJ) const
106 {
107 static_assert(i != j, "Domain i cannot be coupled to itself!");
108 static_assert(AlwaysFalse<Dune::index_constant<i>>::value,
109 "The coupling manager does not implement the couplingStencil() function" );
110
111 return CouplingStencil(); // suppress compiler warning of function having no return statement
112 }
113
114 /*!
115 * \brief returns an iterable container of all indices of degrees of freedom of domain j
116 * that couple with / influence the residual of the given subcontrolvolume face of domain i
117 *
118 * \param domainI the domain index of domain i
119 * \param scvfI the coupled subcontrolvolume face of domain í
120 * \param domainJ the domain index of domain j
121
122 * \note this is a specialization for getting the indices of the coupled cellcentered dofs
123 */
124 358368 const CouplingStencil& couplingStencil(Dune::index_constant<faceIdx> domainI,
125 const SubControlVolumeFace& scvfI,
126 Dune::index_constant<cellCenterIdx> domainJ) const
127 {
128
1/2
✗ Branch 0 not taken.
✓ Branch 1 taken 358368 times.
358368 const auto& connectivityMap = this->problem(domainI).gridGeometry().connectivityMap();
129 358368 return connectivityMap(domainI, domainJ, scvfI.index());
130 }
131
132 /*!
133 * \copydoc Dumux::CouplingManager::evalCouplingResidual()
134 *
135 * \note this is a specialization for calculating the coupled residual for cellcentered dofs
136 * w.r.t. staggered face dof changes
137 */
138 template<class LocalAssemblerI, std::size_t j>
139 34252 decltype(auto) evalCouplingResidual(Dune::index_constant<cellCenterIdx> domainI,
140 const LocalAssemblerI& localAssemblerI,
141 Dune::index_constant<j> domainJ,
142 std::size_t dofIdxGlobalJ) const
143 {
144 static_assert(domainI != domainJ, "Domain i cannot be coupled to itself!");
145
2/4
✓ Branch 0 taken 3620 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 5200 times.
✗ Branch 3 not taken.
37872 return localAssemblerI.evalLocalResidualForCellCenter();
146 }
147
148 /*!
149 * \brief evaluates the face residual of a coupled face of domain i which depends on the variables
150 * at the degree of freedom with index dofIdxGlobalJ of domain j
151 *
152 * \param domainI the domain index of domain i
153 * \param scvfI the subcontrol volume face whose residual shall be evaluated of domain i
154 * \param localAssemblerI the local assembler assembling the element residual of an element of domain i
155 * \param domainJ the domain index of domain j
156 * \param dofIdxGlobalJ the index of the degree of freedom of domain j which has an influence on the element residual of domain i
157 *
158 * \note the default implementation evaluates the complete face residual
159 * if only certain terms of the residual of the residual are coupled
160 * to dof with index dofIdxGlobalJ the function can be overloaded in the coupling manager
161 * \return the face residual
162 */
163 template<class LocalAssemblerI, std::size_t j>
164 34252 decltype(auto) evalCouplingResidual(Dune::index_constant<faceIdx> domainI,
165 const SubControlVolumeFace& scvfI,
166 const LocalAssemblerI& localAssemblerI,
167 Dune::index_constant<j> domainJ,
168 std::size_t dofIdxGlobalJ) const
169 {
170 static_assert(domainI != domainJ, "Domain i cannot be coupled to itself!");
171
2/4
✓ Branch 0 taken 11340 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 22912 times.
✗ Branch 3 not taken.
45592 return localAssemblerI.evalLocalResidualForFace(scvfI);
172 }
173
174 /*!
175 * \brief return the numeric epsilon used for deflecting primary variables of coupled domain i.
176 * \note specialization for non-staggered schemes
177 */
178 template<std::size_t i, typename std::enable_if_t<(GridGeometry<i>::discMethod != DiscretizationMethods::staggered), int> = 0>
179 34 decltype(auto) numericEpsilon(Dune::index_constant<i> id,
180 const std::string& paramGroup) const
181 {
182
2/4
✓ Branch 1 taken 17 times.
✗ Branch 2 not taken.
✓ Branch 5 taken 17 times.
✗ Branch 6 not taken.
34 return ParentType::numericEpsilon(id, paramGroup);
183 }
184
185 /*!
186 * \brief return the numeric epsilon used for deflecting primary variables of coupled domain i.
187 * \note specialization for non-staggered schemes
188 */
189 template<std::size_t i, typename std::enable_if_t<(GridGeometry<i>::discMethod == DiscretizationMethods::staggered), int> = 0>
190 238 decltype(auto) numericEpsilon(Dune::index_constant<i>,
191 const std::string& paramGroup) const
192 {
193 238 constexpr std::size_t numEqCellCenter = Traits::template SubDomain<cellCenterIdx>::PrimaryVariables::dimension;
194 238 constexpr std::size_t numEqFace = Traits::template SubDomain<faceIdx>::PrimaryVariables::dimension;
195 238 constexpr bool isCellCenter = GridGeometry<i>::isCellCenter();
196 238 constexpr std::size_t numEq = isCellCenter ? numEqCellCenter : numEqFace;
197 238 constexpr auto prefix = isCellCenter ? "CellCenter" : "Face";
198
199 try {
200 238 if(paramGroup.empty())
201 194 return NumericEpsilon<typename Traits::Scalar, numEq>(prefix);
202 else
203 228 return NumericEpsilon<typename Traits::Scalar, numEq>(paramGroup + "." + prefix);
204 }
205 catch (Dune::RangeError& e)
206 {
207 DUNE_THROW(Dumux::ParameterException, "For the staggered model, you can to specify \n\n"
208 " CellCenter.Assembly.NumericDifference.PriVarMagnitude = mCC\n"
209 " Face.Assembly.NumericDifference.PriVarMagnitude = mFace\n"
210 " CellCenter.Assembly.NumericDifference.BaseEpsilon = eCC_0 ... eCC_numEqCellCenter-1\n"
211 " Face.Assembly.NumericDifference.BaseEpsilon = eFace_0 ... eFace_numEqFace-1\n\n"
212 "Wrong number of values set for " << prefix << " (has " << numEq << " primary variable(s))\n\n" << e);
213 }
214
215 }
216
217 };
218
219 } //end namespace Dumux
220
221 #endif
222