3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
extendedsourcestencil.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_EMBEDDED_EXTENDEDSOURCESTENCIL_HH
26#define DUMUX_MULTIDOMAIN_EMBEDDED_EXTENDEDSOURCESTENCIL_HH
27
28#include <vector>
29
30#include <dune/common/indices.hh>
31
36
38
44template<class CouplingManager>
46{
47 using MDTraits = typename CouplingManager::MultiDomainTraits;
48 using Scalar = typename MDTraits::Scalar;
49
50 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
51 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
52 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
53 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
54
55 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
56 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
57
58 template<std::size_t id>
59 static constexpr bool isBox()
60 { return GridGeometry<id>::discMethod == DiscretizationMethod::box; }
61public:
68 template<std::size_t id, class JacobianPattern>
69 void extendJacobianPattern(const CouplingManager& couplingManager, Dune::index_constant<id> domainI, JacobianPattern& pattern) const
70 {
71 // add additional dof dependencies
72 for (const auto& element : elements(couplingManager.gridView(domainI)))
73 {
74 const auto& dofs = extendedSourceStencil_(couplingManager, domainI, element);
75 if constexpr (isBox<domainI>())
76 {
77 for (int i = 0; i < element.subEntities(GridView<domainI>::dimension); ++i)
78 for (const auto globalJ : dofs)
79 pattern.add(couplingManager.problem(domainI).gridGeometry().vertexMapper().subIndex(element, i, GridView<domainI>::dimension), globalJ);
80 }
81 else
82 {
83 const auto globalI = couplingManager.problem(domainI).gridGeometry().elementMapper().index(element);
84 for (const auto globalJ : dofs)
85 pattern.add(globalI, globalJ);
86 }
87 }
88 }
89
97 template<std::size_t i, class LocalAssemblerI, class SolutionVector, class JacobianMatrixDiagBlock, class GridVariables>
99 Dune::index_constant<i> domainI,
100 const LocalAssemblerI& localAssemblerI,
101 const SolutionVector& curSol,
102 JacobianMatrixDiagBlock& A,
103 GridVariables& gridVariables) const
104 {
105 constexpr auto numEq = std::decay_t<decltype(curSol[domainI][0])>::size();
106 const auto& elementI = localAssemblerI.element();
107
108 // only do something if we have an extended stencil
109 if (extendedSourceStencil_(couplingManager, domainI, elementI).empty())
110 return;
111
112 // compute the undeflected residual (source only!)
113 const auto origResidual = localAssemblerI.evalLocalSourceResidual(elementI);
114
115 // compute derivate for all additional dofs in the circle stencil
116 for (const auto dofIndex : extendedSourceStencil_(couplingManager, domainI, elementI))
117 {
118 auto partialDerivs = origResidual;
119 const auto origPriVars = curSol[domainI][dofIndex];
120
121 // calculate derivatives w.r.t to the privars at the dof at hand
122 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
123 {
124 // reset partial derivatives
125 partialDerivs = 0.0;
126
127 auto evalResiduals = [&](Scalar priVar)
128 {
129 // update the coupling context (solution vector and recompute element residual)
130 auto priVars = origPriVars;
131 priVars[pvIdx] = priVar;
132 couplingManager.updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, priVars, pvIdx);
133 return localAssemblerI.evalLocalSourceResidual(elementI);
134 };
135
136 // derive the residuals numerically
137 static const int numDiffMethod = getParam<int>("Assembly.NumericDifferenceMethod");
138 NumericDifferentiation::partialDerivative(evalResiduals, curSol[domainI][dofIndex][pvIdx],
139 partialDerivs, origResidual, numDiffMethod);
140
141 // update the global stiffness matrix with the current partial derivatives
142 for (const auto& scvJ : scvs(localAssemblerI.fvGeometry()))
143 {
144 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
145 {
146 // A[i][col][eqIdx][pvIdx] is the rate of change of
147 // the residual of equation 'eqIdx' at dof 'i'
148 // depending on the primary variable 'pvIdx' at dof
149 // 'col'.
150 A[scvJ.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scvJ.indexInElement()][eqIdx];
151 }
152 }
153
154 // restore the original coupling context
155 couplingManager.updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, origPriVars, pvIdx);
156 }
157 }
158 }
159
161 void clear() { sourceStencils_.clear(); }
162
164 typename CouplingManager::template CouplingStencils<bulkIdx>& stencil()
165 { return sourceStencils_; }
166
167private:
169 template<std::size_t id>
170 const auto& extendedSourceStencil_(const CouplingManager& couplingManager, Dune::index_constant<id> domainId, const Element<id>& element) const
171 {
172 if constexpr (domainId == bulkIdx)
173 {
174 const auto bulkElementIdx = couplingManager.problem(bulkIdx).gridGeometry().elementMapper().index(element);
175 if (sourceStencils_.count(bulkElementIdx))
176 return sourceStencils_.at(bulkElementIdx);
177 }
178
179 return couplingManager.emptyStencil(domainId);
180 }
181
183 typename CouplingManager::template CouplingStencils<bulkIdx> sourceStencils_;
184};
185
186} // end namespace Dumux::EmbeddedCoupling
187
188#endif
A class for numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
void updateCouplingContext(Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, Dune::index_constant< j > domainJ, std::size_t dofIdxGlobalJ, const PrimaryVariables< j > &priVarsJ, int pvIdxJ)
updates all data and variables that are necessary to evaluate the residual of the element of domain i...
Definition: multidomain/couplingmanager.hh:152
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
Definition: circlepoints.hh:36
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with repect to a function parameter.
Definition: numericdifferentiation.hh:61
Definition: common/properties.hh:101
Definition: multidomain/couplingmanager.hh:46
const Problem< i > & problem(Dune::index_constant< i > domainIdx) const
Return a reference to the sub problem.
Definition: multidomain/couplingmanager.hh:264
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:46
void extendJacobianPattern(const CouplingManager &couplingManager, Dune::index_constant< id > domainI, JacobianPattern &pattern) const
extend the jacobian pattern of the diagonal block of domain i by those entries that are not already i...
Definition: extendedsourcestencil.hh:69
void evalAdditionalDomainDerivatives(CouplingManager &couplingManager, Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, const SolutionVector &curSol, JacobianMatrixDiagBlock &A, GridVariables &gridVariables) const
evaluate additional derivatives of the element residual of a domain with respect to dofs in the same ...
Definition: extendedsourcestencil.hh:98
void clear()
clear the internal data
Definition: extendedsourcestencil.hh:161
CouplingManager::template CouplingStencils< bulkIdx > & stencil()
return a reference to the stencil
Definition: extendedsourcestencil.hh:164
Declares all properties used in Dumux.