version 3.8
1p/incompressiblelocalresidual.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_1P_INCOMPRESSIBLE_LOCAL_RESIDUAL_HH
15#define DUMUX_1P_INCOMPRESSIBLE_LOCAL_RESIDUAL_HH
16
19
20namespace Dumux {
21
27template<class TypeTag>
29{
34 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
37 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
38 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
39 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
40 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
42 using Element = typename GridView::template Codim<0>::Entity;
44 // first index for the mass balance
45 enum { conti0EqIdx = Indices::conti0EqIdx };
46 enum { pressureIdx = Indices::pressureIdx };
47
48 static constexpr int dim = GridView::dimension;
49 static constexpr int dimWorld = GridView::dimensionworld;
50
51public:
52 using ParentType::ParentType;
53
54 template<class PartialDerivativeMatrix>
55 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
56 const Problem& problem,
57 const Element& element,
58 const FVElementGeometry& fvGeometry,
59 const VolumeVariables& curVolVars,
60 const SubControlVolume& scv) const {}
61
62 template<class PartialDerivativeMatrix>
63 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
64 const Problem& problem,
65 const Element& element,
66 const FVElementGeometry& fvGeometry,
67 const VolumeVariables& curVolVars,
68 const SubControlVolume& scv) const
69 {
70 problem.addSourceDerivatives(partialDerivatives, element, fvGeometry, curVolVars, scv);
71 }
72
74 template<class PartialDerivativeMatrices, class T = TypeTag>
75 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
76 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
77 const Problem& problem,
78 const Element& element,
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& curElemVolVars,
81 const ElementFluxVariablesCache& elemFluxVarsCache,
82 const SubControlVolumeFace& scvf) const
83 {
84 static_assert(!FluidSystem::isCompressible(0),
85 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
86 static_assert(FluidSystem::viscosityIsConstant(0),
87 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
88
89 // we know the "upwind factor" is constant, get inner one here and compute derivatives
90 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
91 / curElemVolVars[scvf.insideScvIdx()].viscosity();
92
93 // for network grids we have to do something else
94 if (dim < dimWorld && scvf.numOutsideScvs() > 1)
95 {
96 const auto insideTi = elemFluxVarsCache[scvf].advectionTij();
97 const auto sumTi = [&]
98 {
99 Scalar sumTi(insideTi);
100 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
101 {
102 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
103 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
104 sumTi += outsideFluxVarsCache.advectionTij();
105 }
106 return sumTi;
107 }();
108
109 // add partial derivatives to the respective given matrices
110 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx]
111 += up*insideTi*(1.0 - insideTi/sumTi);
112 for (unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
113 derivativeMatrices[scvf.outsideScvIdx(i)][conti0EqIdx][pressureIdx]
114 += -up*insideTi*elemFluxVarsCache[fvGeometry.flipScvf(scvf.index(), i)].advectionTij()/sumTi;
115
116 }
117 else
118 {
119 const auto deriv = elemFluxVarsCache[scvf].advectionTij()*up;
120 // add partial derivatives to the respective given matrices
121 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx] += deriv;
122 derivativeMatrices[scvf.outsideScvIdx()][conti0EqIdx][pressureIdx] -= deriv;
123 }
124 }
125
127 template<class PartialDerivativeMatrices, class T = TypeTag>
128 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::ccmpfa, void>
129 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
130 const Problem& problem,
131 const Element& element,
132 const FVElementGeometry& fvGeometry,
133 const ElementVolumeVariables& curElemVolVars,
134 const ElementFluxVariablesCache& elemFluxVarsCache,
135 const SubControlVolumeFace& scvf) const
136 {
137 static_assert(!FluidSystem::isCompressible(0),
138 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
139 static_assert(FluidSystem::viscosityIsConstant(0),
140 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
141
142 // we know the "upwind factor" is constant, get inner one here and compute derivatives
143 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
144 / curElemVolVars[scvf.insideScvIdx()].viscosity();
145
146 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
147 const auto localFaceIdx = fluxVarsCache.ivLocalFaceIndex();
148 const auto usesSecondary = fluxVarsCache.usesSecondaryIv();
149 const auto switchSign = fluxVarsCache.advectionSwitchFluxSign();
150
151 const auto& stencil = fluxVarsCache.advectionStencil();
152 const auto& tij = usesSecondary ? fluxVarsCache.advectionSecondaryDataHandle().T()[localFaceIdx]
153 : fluxVarsCache.advectionPrimaryDataHandle().T()[localFaceIdx];
154
155 // We assume same the tij are order as the stencil up to stencil.size()
156 // any contribution of Dirichlet BCs is assumed to be placed afterwards
157 assert(stencil.size() <= tij.size());
158 for (unsigned int i = 0; i < stencil.size();++i)
159 derivativeMatrices[stencil[i]][conti0EqIdx][pressureIdx] += switchSign ? -tij[i]*up
160 : tij[i]*up;
161 }
162
164 template<class JacobianMatrix, class T = TypeTag>
165 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
166 addFluxDerivatives(JacobianMatrix& A,
167 const Problem& problem,
168 const Element& element,
169 const FVElementGeometry& fvGeometry,
170 const ElementVolumeVariables& curElemVolVars,
171 const ElementFluxVariablesCache& elemFluxVarsCache,
172 const SubControlVolumeFace& scvf) const
173 {
174 static_assert(!FluidSystem::isCompressible(0),
175 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
176 static_assert(FluidSystem::viscosityIsConstant(0),
177 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
178
179 using AdvectionType = GetPropType<T, Properties::AdvectionType>;
180 const auto ti = AdvectionType::calculateTransmissibilities(problem,
181 element,
182 fvGeometry,
183 curElemVolVars,
184 scvf,
185 elemFluxVarsCache[scvf]);
186
187 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
188 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
189
190 // we know the "upwind factor" is constant, get inner one here and compute derivatives
191 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
192 / curElemVolVars[scvf.insideScvIdx()].viscosity();
193 for (const auto& scv : scvs(fvGeometry))
194 {
195 auto d = up*ti[scv.indexInElement()];
196 A[insideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] += d;
197 A[outsideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] -= d;
198 }
199 }
200
202 template<class PartialDerivativeMatrices, class T = TypeTag>
203 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
204 addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
205 const Problem& problem,
206 const Element& element,
207 const FVElementGeometry& fvGeometry,
208 const ElementVolumeVariables& curElemVolVars,
209 const ElementFluxVariablesCache& elemFluxVarsCache,
210 const SubControlVolumeFace& scvf) const
211 {
212 // we know the "upwind factor" is constant, get inner one here
213 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
214 / curElemVolVars[scvf.insideScvIdx()].viscosity();
215 const auto deriv = elemFluxVarsCache[scvf].advectionTij()*up;
216
217 // compute and add partial derivative to the respective given matrices
218 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx] += deriv;
219 }
220
222 template<class PartialDerivativeMatrices, class T = TypeTag>
223 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::ccmpfa, void>
224 addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
225 const Problem& problem,
226 const Element& element,
227 const FVElementGeometry& fvGeometry,
228 const ElementVolumeVariables& curElemVolVars,
229 const ElementFluxVariablesCache& elemFluxVarsCache,
230 const SubControlVolumeFace& scvf) const
231 {
232 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
233 }
234
236 template<class PartialDerivativeMatrices>
237 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
238 const Problem& problem,
239 const Element& element,
240 const FVElementGeometry& fvGeometry,
241 const ElementVolumeVariables& curElemVolVars,
242 const ElementFluxVariablesCache& elemFluxVarsCache,
243 const SubControlVolumeFace& scvf) const
244 {
248 }
249};
250
251} // end namespace Dumux
252
253#endif
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...
Definition: porousmediumflow/immiscible/localresidual.hh:28
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...
Definition: 1p/incompressiblelocalresidual.hh:29
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::box, void > addFluxDerivatives(JacobianMatrix &A, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Flux derivatives for the box scheme.
Definition: 1p/incompressiblelocalresidual.hh:166
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Definition: 1p/incompressiblelocalresidual.hh:55
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Definition: 1p/incompressiblelocalresidual.hh:63
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::cctpfa, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Flux derivatives for the cell-centered tpfa scheme.
Definition: 1p/incompressiblelocalresidual.hh:76
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::ccmpfa, void > addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Dirichlet flux derivatives for the cell-centered mpfa scheme.
Definition: 1p/incompressiblelocalresidual.hh:224
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::cctpfa, void > addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Dirichlet flux derivatives for the cell-centered tpfa scheme.
Definition: 1p/incompressiblelocalresidual.hh:204
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethods::ccmpfa, void > addFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Flux derivatives for the cell-centered mpfa scheme.
Definition: 1p/incompressiblelocalresidual.hh:129
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Robin-type flux derivatives.
Definition: 1p/incompressiblelocalresidual.hh:237
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
constexpr CCMpfa ccmpfa
Definition: method.hh:146
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...