version 3.8
2p/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_2P_INCOMPRESSIBLE_TEST_LOCAL_RESIDUAL_HH
15#define DUMUX_2P_INCOMPRESSIBLE_TEST_LOCAL_RESIDUAL_HH
16
17#include <cmath>
18
21
25
28
29namespace Dumux {
30
36template<class TypeTag>
38{
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
46 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
49 using FVElementGeometry = typename GridGeometry::LocalView;
50 using SubControlVolume = typename GridGeometry::SubControlVolume;
51 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
52 using Extrusion = Extrusion_t<GridGeometry>;
54 using Element = typename GridView::template Codim<0>::Entity;
57
58 static constexpr int numPhases = ModelTraits::numFluidPhases();
59 static constexpr int pressureIdx = ModelTraits::Indices::pressureIdx;
60 static constexpr int saturationIdx = ModelTraits::Indices::saturationIdx;
61 static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx;
62
63public:
64 using ParentType::ParentType;
65
79 template<class PartialDerivativeMatrix>
80 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
81 const Problem& problem,
82 const Element& element,
83 const FVElementGeometry& fvGeometry,
84 const VolumeVariables& curVolVars,
85 const SubControlVolume& scv) const
86 {
87 static_assert(!FluidSystem::isCompressible(0),
88 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
89 static_assert(!FluidSystem::isCompressible(1),
90 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
91 static_assert(ModelTraits::numFluidPhases() == 2,
92 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
93 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
94 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p1-s0 formulation!");
95
96 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*curVolVars.porosity();
97 const auto dt = this->timeLoop().timeStepSize();
98
99 // partial derivative of the phase storage terms w.r.t. S_n
100 partialDerivatives[conti0EqIdx+0][saturationIdx] -= poreVolume*curVolVars.density(0)/dt;
101 partialDerivatives[conti0EqIdx+1][saturationIdx] += poreVolume*curVolVars.density(1)/dt;
102 }
103
114 template<class PartialDerivativeMatrix>
115 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
116 const Problem& problem,
117 const Element& element,
118 const FVElementGeometry& fvGeometry,
119 const VolumeVariables& curVolVars,
120 const SubControlVolume& scv) const
121 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
122
137 template<class PartialDerivativeMatrices, class T = TypeTag>
138 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
139 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
140 const Problem& problem,
141 const Element& element,
142 const FVElementGeometry& fvGeometry,
143 const ElementVolumeVariables& curElemVolVars,
144 const ElementFluxVariablesCache& elemFluxVarsCache,
145 const SubControlVolumeFace& scvf) const
146 {
147 static_assert(!FluidSystem::isCompressible(0),
148 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
149 static_assert(!FluidSystem::isCompressible(1),
150 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
151 static_assert(FluidSystem::viscosityIsConstant(0),
152 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
153 static_assert(FluidSystem::viscosityIsConstant(1),
154 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
155 static_assert(ModelTraits::numFluidPhases() == 2,
156 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
157 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
158 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p1-s0 formulation!");
159
161
162 // evaluate the current wetting phase Darcy flux and resulting upwind weights
163 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
164 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
165 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
166 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
167 const auto outsideWeight_w = 1.0 - insideWeight_w;
168 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
169 const auto outsideWeight_n = 1.0 - insideWeight_n;
170
171 // get references to the two participating vol vars & parameters
172 const auto insideScvIdx = scvf.insideScvIdx();
173 const auto outsideScvIdx = scvf.outsideScvIdx();
174 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
175 const auto& insideScv = fvGeometry.scv(insideScvIdx);
176 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
177 const auto& insideVolVars = curElemVolVars[insideScvIdx];
178 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
179
180 const auto insidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
181 element, insideScv, elementSolution<FVElementGeometry>(insideVolVars.priVars())
182 );
183 const auto outsidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
184 outsideElement, outsideScv, elementSolution<FVElementGeometry>(outsideVolVars.priVars())
185 );
186
187 // get references to the two participating derivative matrices
188 auto& dI_dI = derivativeMatrices[insideScvIdx];
189 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
190
191 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
192 static const auto rho_w = insideVolVars.density(0);
193 static const auto rho_n = insideVolVars.density(1);
194 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
195 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
196 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
197 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
198 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
199 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
200
201 // derivative w.r.t. to Sn is the negative of the one w.r.t. Sw
202 const auto insideSw = insideVolVars.saturation(0);
203 const auto outsideSw = outsideVolVars.saturation(0);
204 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
205 const auto dKrw_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrw_dsw(outsideSw);
206 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
207 const auto dKrn_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrn_dsw(outsideSw);
208 const auto dpc_dSn_inside = -1.0*insidefluidMatrixInteraction.dpc_dsw(insideSw);
209 const auto dpc_dSn_outside = -1.0*outsidefluidMatrixInteraction.dpc_dsw(outsideSw);
210
211 const auto tij = elemFluxVarsCache[scvf].advectionTij();
212
213 // precalculate values
214 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
215 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
216 const auto rho_mu_flux_w = rhow_muw*flux_w;
217 const auto rho_mu_flux_n = rhon_mun*flux_n;
218 const auto tij_up_w = tij*up_w;
219 const auto tij_up_n = tij*up_n;
220
221 // partial derivative of the wetting phase flux w.r.t. p_w
222 dI_dI[conti0EqIdx+0][pressureIdx] += tij_up_w;
223 dI_dJ[conti0EqIdx+0][pressureIdx] -= tij_up_w;
224
225 // partial derivative of the wetting phase flux w.r.t. S_n
226 dI_dI[conti0EqIdx+0][saturationIdx] += rho_mu_flux_w*dKrw_dSn_inside*insideWeight_w;
227 dI_dJ[conti0EqIdx+0][saturationIdx] += rho_mu_flux_w*dKrw_dSn_outside*outsideWeight_w;
228
229 // partial derivative of the nonwetting phase flux w.r.t. p_w
230 dI_dI[conti0EqIdx+1][pressureIdx] += tij_up_n;
231 dI_dJ[conti0EqIdx+1][pressureIdx] -= tij_up_n;
232
233 // partial derivative of the nonwetting phase flux w.r.t. S_n (relative permeability derivative contribution)
234 dI_dI[conti0EqIdx+1][saturationIdx] += rho_mu_flux_n*dKrn_dSn_inside*insideWeight_n;
235 dI_dJ[conti0EqIdx+1][saturationIdx] += rho_mu_flux_n*dKrn_dSn_outside*outsideWeight_n;
236
237 // partial derivative of the nonwetting phase flux w.r.t. S_n (capillary pressure derivative contribution)
238 dI_dI[conti0EqIdx+1][saturationIdx] += tij_up_n*dpc_dSn_inside;
239 dI_dJ[conti0EqIdx+1][saturationIdx] -= tij_up_n*dpc_dSn_outside;
240 }
241
256 template<class JacobianMatrix, class T = TypeTag>
257 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
258 addFluxDerivatives(JacobianMatrix& A,
259 const Problem& problem,
260 const Element& element,
261 const FVElementGeometry& fvGeometry,
262 const ElementVolumeVariables& curElemVolVars,
263 const ElementFluxVariablesCache& elemFluxVarsCache,
264 const SubControlVolumeFace& scvf) const
265 {
266 static_assert(!FluidSystem::isCompressible(0),
267 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
268 static_assert(!FluidSystem::isCompressible(1),
269 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
270 static_assert(FluidSystem::viscosityIsConstant(0),
271 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
272 static_assert(FluidSystem::viscosityIsConstant(1),
273 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
274 static_assert(ModelTraits::numFluidPhases() == 2,
275 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
276 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
277 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p0-s1 formulation!");
278
280
281 // evaluate the current wetting phase Darcy flux and resulting upwind weights
282 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
283 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
284 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
285 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
286 const auto outsideWeight_w = 1.0 - insideWeight_w;
287 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
288 const auto outsideWeight_n = 1.0 - insideWeight_n;
289
290 // get references to the two participating vol vars & parameters
291 const auto insideScvIdx = scvf.insideScvIdx();
292 const auto outsideScvIdx = scvf.outsideScvIdx();
293 const auto& insideScv = fvGeometry.scv(insideScvIdx);
294 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
295 const auto& insideVolVars = curElemVolVars[insideScv];
296 const auto& outsideVolVars = curElemVolVars[outsideScv];
297
298 const auto elemSol = elementSolution(element, curElemVolVars, fvGeometry);
299
300 const auto insidefluidMatrixInteraction
301 = problem.spatialParams().fluidMatrixInteraction(element, insideScv, elemSol);
302 const auto outsidefluidMatrixInteraction
303 = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, elemSol);
304
305 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
306 static const auto rho_w = insideVolVars.density(0);
307 static const auto rho_n = insideVolVars.density(1);
308 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
309 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
310 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
311 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
312 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
313 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
314
315 // let the Law for the advective fluxes calculate the transmissibilities
316 const auto ti = AdvectionType::calculateTransmissibilities(problem,
317 element,
318 fvGeometry,
319 curElemVolVars,
320 scvf,
321 elemFluxVarsCache[scvf]);
322
323 // get the rows of the jacobian matrix for the inside/outside scv
324 auto& dI_dJ_inside = A[insideScv.dofIndex()];
325 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
326
327 // precalculate values
328 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
329 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
330 const auto rho_mu_flux_w = rhow_muw*flux_w;
331 const auto rho_mu_flux_n = rhon_mun*flux_n;
332
333 // add the partial derivatives w.r.t all scvs in the element
334 for (const auto& scvJ : scvs(fvGeometry))
335 {
336 const auto globalJ = scvJ.dofIndex();
337 const auto localJ = scvJ.indexInElement();
338
339 // the transmissibily associated with the scvJ
340 const auto tj = ti[localJ];
341
342 // partial derivative of the wetting phase flux w.r.t. p_w
343 const auto tj_up_w = tj*up_w;
344 dI_dJ_inside[globalJ][conti0EqIdx+0][pressureIdx] += tj_up_w;
345 dI_dJ_outside[globalJ][conti0EqIdx+0][pressureIdx] -= tj_up_w;
346
347 // partial derivative of the nonwetting phase flux w.r.t. p_w
348 const auto tj_up_n = tj*up_n;
349 dI_dJ_inside[globalJ][conti0EqIdx+1][pressureIdx] += tj_up_n;
350 dI_dJ_outside[globalJ][conti0EqIdx+1][pressureIdx] -= tj_up_n;
351
352 // partial derivatives w.r.t. S_n (are the negative of those w.r.t sw)
353 // relative permeability contributions only for inside/outside
354 if (localJ == insideScvIdx)
355 {
356 // partial derivative of the wetting phase flux w.r.t. S_n
357 const auto insideSw = insideVolVars.saturation(0);
358 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
359 const auto dFluxW_dSnJ = rho_mu_flux_w*dKrw_dSn_inside*insideWeight_w;
360 dI_dJ_inside[globalJ][conti0EqIdx+0][saturationIdx] += dFluxW_dSnJ;
361 dI_dJ_outside[globalJ][conti0EqIdx+0][saturationIdx] -= dFluxW_dSnJ;
362
363 // partial derivative of the nonwetting phase flux w.r.t. S_n (k_rn contribution)
364 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
365 const auto dFluxN_dSnJ_krn = rho_mu_flux_n*dKrn_dSn_inside*insideWeight_n;
366 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_krn;
367 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_krn;
368
369 // partial derivative of the nonwetting phase flux w.r.t. S_n (p_c contribution)
370 const auto dFluxN_dSnJ_pc = -1.0*tj_up_n*insidefluidMatrixInteraction.dpc_dsw(insideSw);
371 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
372 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
373 }
374 else if (localJ == outsideScvIdx)
375 {
376 // see comments for (globalJ == insideScvIdx)
377 const auto outsideSw = outsideVolVars.saturation(0);
378 const auto dKrw_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrw_dsw(outsideSw);
379 const auto dFluxW_dSnJ = rho_mu_flux_w*dKrw_dSn_outside*outsideWeight_w;
380 dI_dJ_inside[globalJ][conti0EqIdx+0][saturationIdx] += dFluxW_dSnJ;
381 dI_dJ_outside[globalJ][conti0EqIdx+0][saturationIdx] -= dFluxW_dSnJ;
382
383 const auto dKrn_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrn_dsw(outsideSw);
384 const auto dFluxN_dSnJ_krn = rho_mu_flux_n*dKrn_dSn_outside*outsideWeight_n;
385 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_krn;
386 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_krn;
387
388 const auto dFluxN_dSnJ_pc = -1.0*tj_up_n*outsidefluidMatrixInteraction.dpc_dsw(outsideSw);
389 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
390 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
391 }
392 else
393 {
394 const auto& fluidMatrixInteractionJ
395 = problem.spatialParams().fluidMatrixInteraction(element, scvJ, elemSol);
396
397 const auto swJ = curElemVolVars[scvJ].saturation(0);
398 const auto dFluxN_dSnJ_pc = tj_up_n*fluidMatrixInteractionJ.dpc_dsw(swJ);
399 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
400 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
401 }
402 }
403 }
404
419 template<class PartialDerivativeMatrices>
420 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
421 const Problem& problem,
422 const Element& element,
423 const FVElementGeometry& fvGeometry,
424 const ElementVolumeVariables& curElemVolVars,
425 const ElementFluxVariablesCache& elemFluxVarsCache,
426 const SubControlVolumeFace& scvf) const
427 {
429
430 // evaluate the current wetting phase Darcy flux and resulting upwind weights
431 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
432 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
433 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
434 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
435 const auto outsideWeight_w = 1.0 - insideWeight_w;
436 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
437 const auto outsideWeight_n = 1.0 - insideWeight_n;
438
439 // get references to the two participating vol vars & parameters
440 const auto insideScvIdx = scvf.insideScvIdx();
441 const auto& insideScv = fvGeometry.scv(insideScvIdx);
442 const auto& insideVolVars = curElemVolVars[insideScvIdx];
443 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
444
445 const auto insidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
446 element, insideScv, elementSolution<FVElementGeometry>(insideVolVars.priVars())
447 );
448
449 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
450 static const auto rho_w = insideVolVars.density(0);
451 static const auto rho_n = insideVolVars.density(1);
452 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
453 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
454 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
455 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
456 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
457 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
458
459 // get reference to the inside derivative matrix
460 auto& dI_dI = derivativeMatrices[insideScvIdx];
461
462 // derivative w.r.t. to Sn is the negative of the one w.r.t. Sw
463 const auto insideSw = insideVolVars.saturation(0);
464 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
465 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
466 const auto dpc_dSn_inside = -1.0*insidefluidMatrixInteraction.dpc_dsw(insideSw);
467
468 const auto tij = elemFluxVarsCache[scvf].advectionTij();
469 // partial derivative of the wetting phase flux w.r.t. p_w
470 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
471 dI_dI[conti0EqIdx+0][pressureIdx] += tij*up_w;
472
473 // partial derivative of the wetting phase flux w.r.t. S_n
474 dI_dI[conti0EqIdx+0][saturationIdx] += rhow_muw*flux_w*dKrw_dSn_inside*insideWeight_w;
475
476 // partial derivative of the nonwetting phase flux w.r.t. p_w
477 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
478 dI_dI[conti0EqIdx+1][pressureIdx] += tij*up_n;
479
480 // partial derivative of the nonwetting phase flux w.r.t. S_n (relative permeability derivative contribution)
481 dI_dI[conti0EqIdx+1][saturationIdx] += rhon_mun*flux_n*dKrn_dSn_inside*insideWeight_n;
482
483 // partial derivative of the nonwetting phase flux w.r.t. S_n (capillary pressure derivative contribution)
484 dI_dI[conti0EqIdx+1][saturationIdx] += tij*dpc_dSn_inside*up_n;
485 }
486
498 template<class PartialDerivativeMatrices>
499 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
500 const Problem& problem,
501 const Element& element,
502 const FVElementGeometry& fvGeometry,
503 const ElementVolumeVariables& curElemVolVars,
504 const ElementFluxVariablesCache& elemFluxVarsCache,
505 const SubControlVolumeFace& scvf) const
506 { /* TODO maybe forward to problem for the user to implement the Robin derivatives?*/ }
507};
508
509} // end namespace Dumux
510
511#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 two-phase, incompressible test pro...
Definition: 2p/incompressiblelocalresidual.hh:38
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
Adds flux derivatives for wetting and nonwetting phase for cell-centered FVM using TPFA.
Definition: 2p/incompressiblelocalresidual.hh:139
void addStorageDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Adds storage derivatives for wetting and nonwetting phase.
Definition: 2p/incompressiblelocalresidual.hh:80
void addRobinFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Adds Robin flux derivatives for wetting and nonwetting phase.
Definition: 2p/incompressiblelocalresidual.hh:499
void addSourceDerivatives(PartialDerivativeMatrix &partialDerivatives, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const VolumeVariables &curVolVars, const SubControlVolume &scv) const
Adds source derivatives for wetting and nonwetting phase.
Definition: 2p/incompressiblelocalresidual.hh:115
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
Adds flux derivatives for wetting and nonwetting phase for box method.
Definition: 2p/incompressiblelocalresidual.hh:258
void addCCDirichletFluxDerivatives(PartialDerivativeMatrices &derivativeMatrices, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Adds cell-centered Dirichlet flux derivatives for wetting and nonwetting phase.
Definition: 2p/incompressiblelocalresidual.hh:420
Defines all properties used in Dumux.
Element solution classes and factory functions.
Helper classes to compute the integration elements.
Defines an enumeration for the formulations accepted by the two-phase model.
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::cctpfa||GridGeometry::discMethod==DiscretizationMethods::ccmpfa, CCElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for cell-centered schemes.
Definition: cellcentered/elementsolution.hh:101
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
@ p0s1
first phase pressure and second phase saturation as primary variables
The available discretization methods in Dumux.
constexpr CCTpfa cctpfa
Definition: method.hh:145
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...