26#ifndef DUMUX_1P_INCOMPRESSIBLE_LOCAL_RESIDUAL_HH
27#define DUMUX_1P_INCOMPRESSIBLE_LOCAL_RESIDUAL_HH
39template<
class TypeTag>
51 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
54 using Element =
typename GridView::template Codim<0>::Entity;
57 enum { conti0EqIdx = Indices::conti0EqIdx };
58 enum { pressureIdx = Indices::pressureIdx };
60 static constexpr int dim = GridView::dimension;
61 static constexpr int dimWorld = GridView::dimensionworld;
64 using ParentType::ParentType;
66 template<
class PartialDerivativeMatrix>
68 const Problem& problem,
69 const Element& element,
70 const FVElementGeometry& fvGeometry,
71 const VolumeVariables& curVolVars,
72 const SubControlVolume& scv)
const {}
74 template<
class PartialDerivativeMatrix>
76 const Problem& problem,
77 const Element& element,
78 const FVElementGeometry& fvGeometry,
79 const VolumeVariables& curVolVars,
80 const SubControlVolume& scv)
const
82 problem.addSourceDerivatives(partialDerivatives, element, fvGeometry, curVolVars, scv);
86 template<
class PartialDerivativeMatrices,
class T = TypeTag>
89 const Problem& problem,
90 const Element& element,
91 const FVElementGeometry& fvGeometry,
92 const ElementVolumeVariables& curElemVolVars,
93 const ElementFluxVariablesCache& elemFluxVarsCache,
94 const SubControlVolumeFace& scvf)
const
96 static_assert(!FluidSystem::isCompressible(0),
97 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
98 static_assert(FluidSystem::viscosityIsConstant(0),
99 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
102 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
103 / curElemVolVars[scvf.insideScvIdx()].viscosity();
106 if (dim < dimWorld && scvf.numOutsideScvs() > 1)
108 const auto insideTi = elemFluxVarsCache[scvf].advectionTij();
109 const auto sumTi = [&]
111 Scalar sumTi(insideTi);
112 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
114 const auto& flippedScvf = fvGeometry.flipScvf(scvf.index(), i);
115 const auto& outsideFluxVarsCache = elemFluxVarsCache[flippedScvf];
116 sumTi += outsideFluxVarsCache.advectionTij();
122 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx]
123 += up*insideTi*(1.0 - insideTi/sumTi);
124 for (
unsigned int i = 0; i < scvf.numOutsideScvs(); ++i)
125 derivativeMatrices[scvf.outsideScvIdx(i)][conti0EqIdx][pressureIdx]
126 += -up*insideTi*elemFluxVarsCache[fvGeometry.flipScvf(scvf.index(), i)].advectionTij()/sumTi;
131 const auto deriv = elemFluxVarsCache[scvf].advectionTij()*up;
133 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx] += deriv;
134 derivativeMatrices[scvf.outsideScvIdx()][conti0EqIdx][pressureIdx] -= deriv;
139 template<
class PartialDerivativeMatrices,
class T = TypeTag>
142 const Problem& problem,
143 const Element& element,
144 const FVElementGeometry& fvGeometry,
145 const ElementVolumeVariables& curElemVolVars,
146 const ElementFluxVariablesCache& elemFluxVarsCache,
147 const SubControlVolumeFace& scvf)
const
149 static_assert(!FluidSystem::isCompressible(0),
150 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
151 static_assert(FluidSystem::viscosityIsConstant(0),
152 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
155 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
156 / curElemVolVars[scvf.insideScvIdx()].viscosity();
158 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
159 const auto localFaceIdx = fluxVarsCache.ivLocalFaceIndex();
160 const auto usesSecondary = fluxVarsCache.usesSecondaryIv();
161 const auto switchSign = fluxVarsCache.advectionSwitchFluxSign();
163 const auto& stencil = fluxVarsCache.advectionStencil();
164 const auto& tij = usesSecondary ? fluxVarsCache.advectionSecondaryDataHandle().T()[localFaceIdx]
165 : fluxVarsCache.advectionPrimaryDataHandle().T()[localFaceIdx];
169 assert(stencil.size() <= tij.size());
170 for (
unsigned int i = 0; i < stencil.size();++i)
171 derivativeMatrices[stencil[i]][conti0EqIdx][pressureIdx] += switchSign ? -tij[i]*up
176 template<
class JacobianMatrix,
class T = TypeTag>
179 const Problem& problem,
180 const Element& element,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& curElemVolVars,
183 const ElementFluxVariablesCache& elemFluxVarsCache,
184 const SubControlVolumeFace& scvf)
const
186 static_assert(!FluidSystem::isCompressible(0),
187 "1p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
188 static_assert(FluidSystem::viscosityIsConstant(0),
189 "1p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
192 const auto ti = AdvectionType::calculateTransmissibilities(problem,
197 elemFluxVarsCache[scvf]);
199 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
200 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
203 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
204 / curElemVolVars[scvf.insideScvIdx()].viscosity();
205 for (
const auto& scv : scvs(fvGeometry))
207 auto d = up*ti[scv.indexInElement()];
208 A[insideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] += d;
209 A[outsideScv.dofIndex()][scv.dofIndex()][conti0EqIdx][pressureIdx] -= d;
214 template<
class PartialDerivativeMatrices,
class T = TypeTag>
217 const Problem& problem,
218 const Element& element,
219 const FVElementGeometry& fvGeometry,
220 const ElementVolumeVariables& curElemVolVars,
221 const ElementFluxVariablesCache& elemFluxVarsCache,
222 const SubControlVolumeFace& scvf)
const
225 static const Scalar up = curElemVolVars[scvf.insideScvIdx()].density()
226 / curElemVolVars[scvf.insideScvIdx()].viscosity();
227 const auto deriv = elemFluxVarsCache[scvf].advectionTij()*up;
230 derivativeMatrices[scvf.insideScvIdx()][conti0EqIdx][pressureIdx] += deriv;
234 template<
class PartialDerivativeMatrices,
class T = TypeTag>
237 const Problem& problem,
238 const Element& element,
239 const FVElementGeometry& fvGeometry,
240 const ElementVolumeVariables& curElemVolVars,
241 const ElementFluxVariablesCache& elemFluxVarsCache,
242 const SubControlVolumeFace& scvf)
const
244 addFluxDerivatives(derivativeMatrices, problem, element, fvGeometry, curElemVolVars, elemFluxVarsCache, scvf);
248 template<
class PartialDerivativeMatrices>
250 const Problem& problem,
251 const Element& element,
252 const FVElementGeometry& fvGeometry,
253 const ElementVolumeVariables& curElemVolVars,
254 const ElementFluxVariablesCache& elemFluxVarsCache,
255 const SubControlVolumeFace& scvf)
const
The available discretization methods in Dumux.
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
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...
Definition: 1p/incompressiblelocalresidual.hh:41
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:141
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:216
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Definition: 1p/incompressiblelocalresidual.hh:67
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:236
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:88
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Definition: 1p/incompressiblelocalresidual.hh:75
std::enable_if_t< GetPropType< T, Properties::GridGeometry >::discMethod==DiscretizationMethod::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:178
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:249
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...
Definition: porousmediumflow/immiscible/localresidual.hh:39
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...