3.6-git
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
40
46template<class CouplingManager>
48{
49 using MDTraits = typename CouplingManager::MultiDomainTraits;
50 using Scalar = typename MDTraits::Scalar;
51
52 template<std::size_t id> using SubDomainTypeTag = typename MDTraits::template SubDomain<id>::TypeTag;
53 template<std::size_t id> using GridGeometry = GetPropType<SubDomainTypeTag<id>, Properties::GridGeometry>;
54 template<std::size_t id> using GridView = typename GridGeometry<id>::GridView;
55 template<std::size_t id> using Element = typename GridView<id>::template Codim<0>::Entity;
56
57 static constexpr auto bulkIdx = typename MDTraits::template SubDomain<0>::Index();
58 static constexpr auto lowDimIdx = typename MDTraits::template SubDomain<1>::Index();
59
60 template<std::size_t id>
61 static constexpr bool isBox()
62 { return GridGeometry<id>::discMethod == DiscretizationMethods::box; }
63public:
70 template<std::size_t id, class JacobianPattern>
71 void extendJacobianPattern(const CouplingManager& couplingManager, Dune::index_constant<id> domainI, JacobianPattern& pattern) const
72 {
73 // add additional dof dependencies
74 for (const auto& element : elements(couplingManager.gridView(domainI)))
75 {
76 const auto& dofs = extendedSourceStencil_(couplingManager, domainI, element);
77 if constexpr (isBox<domainI>())
78 {
79 for (int i = 0; i < element.subEntities(GridView<domainI>::dimension); ++i)
80 for (const auto globalJ : dofs)
81 pattern.add(couplingManager.problem(domainI).gridGeometry().vertexMapper().subIndex(element, i, GridView<domainI>::dimension), globalJ);
82 }
83 else
84 {
85 const auto globalI = couplingManager.problem(domainI).gridGeometry().elementMapper().index(element);
86 for (const auto globalJ : dofs)
87 pattern.add(globalI, globalJ);
88 }
89 }
90 }
91
99 template<std::size_t i, class LocalAssemblerI, class JacobianMatrixDiagBlock, class GridVariables>
101 Dune::index_constant<i> domainI,
102 const LocalAssemblerI& localAssemblerI,
103 JacobianMatrixDiagBlock& A,
104 GridVariables& gridVariables) const
105 {
106 const auto& curSolI = couplingManager.curSol(domainI);
107 constexpr auto numEq = std::decay_t<decltype(curSolI[0])>::size();
108 const auto& elementI = localAssemblerI.element();
109
110 // only do something if we have an extended stencil
111 if (extendedSourceStencil_(couplingManager, domainI, elementI).empty())
112 return;
113
114 // compute the undeflected residual (source only!)
115 const auto origResidual = localAssemblerI.evalLocalSourceResidual(elementI);
116
117 // compute derivate for all additional dofs in the circle stencil
118 for (const auto dofIndex : extendedSourceStencil_(couplingManager, domainI, elementI))
119 {
120 auto partialDerivs = origResidual;
121 const auto origPriVars = curSolI[dofIndex];
122 auto priVars = origPriVars;
123
124 // calculate derivatives w.r.t to the privars at the dof at hand
125 for (int pvIdx = 0; pvIdx < numEq; pvIdx++)
126 {
127 // reset partial derivatives
128 partialDerivs = 0.0;
129
130 const auto evalResiduals = [&](const Scalar priVar)
131 {
132 // update the coupling context (solution vector and recompute element residual)
133 priVars[pvIdx] = priVar;
134 couplingManager.updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, priVars, pvIdx);
135 return localAssemblerI.evalLocalSourceResidual(elementI);
136 };
137
138 // derive the residuals numerically
139 static const NumericEpsilon<Scalar, numEq> eps_{localAssemblerI.problem().paramGroup()};
140 static const int numDiffMethod = getParamFromGroup<int>(localAssemblerI.problem().paramGroup(), "Assembly.NumericDifferenceMethod");
142 evalResiduals, priVars[pvIdx], partialDerivs, origResidual, eps_(priVars[pvIdx], pvIdx), numDiffMethod
143 );
144
145 // update the global stiffness matrix with the current partial derivatives
146 for (const auto& scvJ : scvs(localAssemblerI.fvGeometry()))
147 {
148 for (int eqIdx = 0; eqIdx < numEq; eqIdx++)
149 {
150 // A[i][col][eqIdx][pvIdx] is the rate of change of
151 // the residual of equation 'eqIdx' at dof 'i'
152 // depending on the primary variable 'pvIdx' at dof
153 // 'col'.
154 A[scvJ.dofIndex()][dofIndex][eqIdx][pvIdx] += partialDerivs[scvJ.indexInElement()][eqIdx];
155 }
156 }
157
158 // restore the current element solution
159 priVars[pvIdx] = origPriVars[pvIdx];
160
161 // restore the original coupling context
162 couplingManager.updateCouplingContext(domainI, localAssemblerI, domainI, dofIndex, origPriVars, pvIdx);
163 }
164 }
165 }
166
168 void clear() { sourceStencils_.clear(); }
169
171 typename CouplingManager::template CouplingStencils<bulkIdx>& stencil()
172 { return sourceStencils_; }
173
175 const typename CouplingManager::template CouplingStencils<bulkIdx>& stencil() const
176 { return sourceStencils_; }
177
178private:
180 template<std::size_t id>
181 const auto& extendedSourceStencil_(const CouplingManager& couplingManager, Dune::index_constant<id> domainId, const Element<id>& element) const
182 {
183 if constexpr (domainId == bulkIdx)
184 {
185 const auto bulkElementIdx = couplingManager.problem(bulkIdx).gridGeometry().elementMapper().index(element);
186 if (auto stencil = sourceStencils_.find(bulkElementIdx); stencil != sourceStencils_.end())
187 return stencil->second;
188 }
189
190 return couplingManager.emptyStencil(domainId);
191 }
192
194 typename CouplingManager::template CouplingStencils<bulkIdx> sourceStencils_;
195};
196
197} // end namespace Dumux::EmbeddedCoupling
198
199#endif
An adapter class for local assemblers using numeric differentiation.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A class for numeric differentiation.
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:195
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Definition: circlepoints.hh:36
A helper class for local assemblers using numeric differentiation to determine the epsilon.
Definition: numericepsilon.hh:41
static void partialDerivative(const Function &function, Scalar x0, FunctionEvalType &derivative, const FunctionEvalType &fx0, const int numericDifferenceMethod=1)
Computes the derivative of a function with respect to a function parameter.
Definition: numericdifferentiation.hh:61
Definition: common/properties.hh:100
The interface of the coupling manager for multi domain problems.
Definition: multidomain/couplingmanager.hh:60
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
A class managing an extended source stencil.
Definition: extendedsourcestencil.hh:48
void evalAdditionalDomainDerivatives(CouplingManager &couplingManager, Dune::index_constant< i > domainI, const LocalAssemblerI &localAssemblerI, 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:100
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:71
const CouplingManager::template CouplingStencils< bulkIdx > & stencil() const
return a const reference to the stencil
Definition: extendedsourcestencil.hh:175
void clear()
clear the internal data
Definition: extendedsourcestencil.hh:168
CouplingManager::template CouplingStencils< bulkIdx > & stencil()
return a reference to the stencil
Definition: extendedsourcestencil.hh:171
Declares all properties used in Dumux.