3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 *****************************************************************************/
26#ifndef DUMUX_2P_INCOMPRESSIBLE_TEST_LOCAL_RESIDUAL_HH
27#define DUMUX_2P_INCOMPRESSIBLE_TEST_LOCAL_RESIDUAL_HH
28
29#include <cmath>
30
33
37
40
41namespace Dumux {
42
48template<class TypeTag>
50{
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
61 using FVElementGeometry = typename GridGeometry::LocalView;
62 using SubControlVolume = typename GridGeometry::SubControlVolume;
63 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
64 using Extrusion = Extrusion_t<GridGeometry>;
66 using Element = typename GridView::template Codim<0>::Entity;
69
70 static constexpr int numPhases = ModelTraits::numFluidPhases();
71 static constexpr int pressureIdx = ModelTraits::Indices::pressureIdx;
72 static constexpr int saturationIdx = ModelTraits::Indices::saturationIdx;
73 static constexpr int conti0EqIdx = ModelTraits::Indices::conti0EqIdx;
74
75public:
76 using ParentType::ParentType;
77
91 template<class PartialDerivativeMatrix>
92 void addStorageDerivatives(PartialDerivativeMatrix& partialDerivatives,
93 const Problem& problem,
94 const Element& element,
95 const FVElementGeometry& fvGeometry,
96 const VolumeVariables& curVolVars,
97 const SubControlVolume& scv) const
98 {
99 static_assert(!FluidSystem::isCompressible(0),
100 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
101 static_assert(!FluidSystem::isCompressible(1),
102 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
103 static_assert(ModelTraits::numFluidPhases() == 2,
104 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
105 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
106 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p1-s0 formulation!");
107
108 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*curVolVars.porosity();
109 const auto dt = this->timeLoop().timeStepSize();
110
111 // partial derivative of the phase storage terms w.r.t. S_n
112 partialDerivatives[conti0EqIdx+0][saturationIdx] -= poreVolume*curVolVars.density(0)/dt;
113 partialDerivatives[conti0EqIdx+1][saturationIdx] += poreVolume*curVolVars.density(1)/dt;
114 }
115
126 template<class PartialDerivativeMatrix>
127 void addSourceDerivatives(PartialDerivativeMatrix& partialDerivatives,
128 const Problem& problem,
129 const Element& element,
130 const FVElementGeometry& fvGeometry,
131 const VolumeVariables& curVolVars,
132 const SubControlVolume& scv) const
133 { /* TODO maybe forward to problem for the user to implement the source derivatives?*/ }
134
149 template<class PartialDerivativeMatrices, class T = TypeTag>
150 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::cctpfa, void>
151 addFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
152 const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& curElemVolVars,
156 const ElementFluxVariablesCache& elemFluxVarsCache,
157 const SubControlVolumeFace& scvf) const
158 {
159 static_assert(!FluidSystem::isCompressible(0),
160 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
161 static_assert(!FluidSystem::isCompressible(1),
162 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
163 static_assert(FluidSystem::viscosityIsConstant(0),
164 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
165 static_assert(FluidSystem::viscosityIsConstant(1),
166 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
167 static_assert(ModelTraits::numFluidPhases() == 2,
168 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
169 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
170 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p1-s0 formulation!");
171
173
174 // evaluate the current wetting phase Darcy flux and resulting upwind weights
175 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
176 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
177 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
178 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
179 const auto outsideWeight_w = 1.0 - insideWeight_w;
180 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
181 const auto outsideWeight_n = 1.0 - insideWeight_n;
182
183 // get references to the two participating vol vars & parameters
184 const auto insideScvIdx = scvf.insideScvIdx();
185 const auto outsideScvIdx = scvf.outsideScvIdx();
186 const auto outsideElement = fvGeometry.gridGeometry().element(outsideScvIdx);
187 const auto& insideScv = fvGeometry.scv(insideScvIdx);
188 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
189 const auto& insideVolVars = curElemVolVars[insideScvIdx];
190 const auto& outsideVolVars = curElemVolVars[outsideScvIdx];
191
192 const auto insidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
193 element, insideScv, elementSolution<FVElementGeometry>(insideVolVars.priVars())
194 );
195 const auto outsidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
196 outsideElement, outsideScv, elementSolution<FVElementGeometry>(outsideVolVars.priVars())
197 );
198
199 // get references to the two participating derivative matrices
200 auto& dI_dI = derivativeMatrices[insideScvIdx];
201 auto& dI_dJ = derivativeMatrices[outsideScvIdx];
202
203 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
204 static const auto rho_w = insideVolVars.density(0);
205 static const auto rho_n = insideVolVars.density(1);
206 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
207 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
208 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
209 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
210 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
211 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
212
213 // derivative w.r.t. to Sn is the negative of the one w.r.t. Sw
214 const auto insideSw = insideVolVars.saturation(0);
215 const auto outsideSw = outsideVolVars.saturation(0);
216 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
217 const auto dKrw_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrw_dsw(outsideSw);
218 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
219 const auto dKrn_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrn_dsw(outsideSw);
220 const auto dpc_dSn_inside = -1.0*insidefluidMatrixInteraction.dpc_dsw(insideSw);
221 const auto dpc_dSn_outside = -1.0*outsidefluidMatrixInteraction.dpc_dsw(outsideSw);
222
223 const auto tij = elemFluxVarsCache[scvf].advectionTij();
224
225 // precalculate values
226 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
227 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
228 const auto rho_mu_flux_w = rhow_muw*flux_w;
229 const auto rho_mu_flux_n = rhon_mun*flux_n;
230 const auto tij_up_w = tij*up_w;
231 const auto tij_up_n = tij*up_n;
232
233 // partial derivative of the wetting phase flux w.r.t. p_w
234 dI_dI[conti0EqIdx+0][pressureIdx] += tij_up_w;
235 dI_dJ[conti0EqIdx+0][pressureIdx] -= tij_up_w;
236
237 // partial derivative of the wetting phase flux w.r.t. S_n
238 dI_dI[conti0EqIdx+0][saturationIdx] += rho_mu_flux_w*dKrw_dSn_inside*insideWeight_w;
239 dI_dJ[conti0EqIdx+0][saturationIdx] += rho_mu_flux_w*dKrw_dSn_outside*outsideWeight_w;
240
241 // partial derivative of the nonwetting phase flux w.r.t. p_w
242 dI_dI[conti0EqIdx+1][pressureIdx] += tij_up_n;
243 dI_dJ[conti0EqIdx+1][pressureIdx] -= tij_up_n;
244
245 // partial derivative of the nonwetting phase flux w.r.t. S_n (relative permeability derivative contribution)
246 dI_dI[conti0EqIdx+1][saturationIdx] += rho_mu_flux_n*dKrn_dSn_inside*insideWeight_n;
247 dI_dJ[conti0EqIdx+1][saturationIdx] += rho_mu_flux_n*dKrn_dSn_outside*outsideWeight_n;
248
249 // partial derivative of the nonwetting phase flux w.r.t. S_n (capillary pressure derivative contribution)
250 dI_dI[conti0EqIdx+1][saturationIdx] += tij_up_n*dpc_dSn_inside;
251 dI_dJ[conti0EqIdx+1][saturationIdx] -= tij_up_n*dpc_dSn_outside;
252 }
253
268 template<class JacobianMatrix, class T = TypeTag>
269 std::enable_if_t<GetPropType<T, Properties::GridGeometry>::discMethod == DiscretizationMethods::box, void>
270 addFluxDerivatives(JacobianMatrix& A,
271 const Problem& problem,
272 const Element& element,
273 const FVElementGeometry& fvGeometry,
274 const ElementVolumeVariables& curElemVolVars,
275 const ElementFluxVariablesCache& elemFluxVarsCache,
276 const SubControlVolumeFace& scvf) const
277 {
278 static_assert(!FluidSystem::isCompressible(0),
279 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
280 static_assert(!FluidSystem::isCompressible(1),
281 "2p/incompressiblelocalresidual.hh: Only incompressible fluids are allowed!");
282 static_assert(FluidSystem::viscosityIsConstant(0),
283 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
284 static_assert(FluidSystem::viscosityIsConstant(1),
285 "2p/incompressiblelocalresidual.hh: Only fluids with constant viscosities are allowed!");
286 static_assert(ModelTraits::numFluidPhases() == 2,
287 "2p/incompressiblelocalresidual.hh: Only two-phase models are allowed!");
288 static_assert(ModelTraits::priVarFormulation() == TwoPFormulation::p0s1,
289 "2p/incompressiblelocalresidual.hh: Analytic differentiation has to be checked for p0-s1 formulation!");
290
292
293 // evaluate the current wetting phase Darcy flux and resulting upwind weights
294 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
295 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
296 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
297 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
298 const auto outsideWeight_w = 1.0 - insideWeight_w;
299 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
300 const auto outsideWeight_n = 1.0 - insideWeight_n;
301
302 // get references to the two participating vol vars & parameters
303 const auto insideScvIdx = scvf.insideScvIdx();
304 const auto outsideScvIdx = scvf.outsideScvIdx();
305 const auto& insideScv = fvGeometry.scv(insideScvIdx);
306 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
307 const auto& insideVolVars = curElemVolVars[insideScv];
308 const auto& outsideVolVars = curElemVolVars[outsideScv];
309
310 const auto elemSol = elementSolution(element, curElemVolVars, fvGeometry);
311
312 const auto insidefluidMatrixInteraction
313 = problem.spatialParams().fluidMatrixInteraction(element, insideScv, elemSol);
314 const auto outsidefluidMatrixInteraction
315 = problem.spatialParams().fluidMatrixInteraction(element, outsideScv, elemSol);
316
317 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
318 static const auto rho_w = insideVolVars.density(0);
319 static const auto rho_n = insideVolVars.density(1);
320 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
321 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
322 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
323 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
324 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
325 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
326
327 // let the Law for the advective fluxes calculate the transmissibilities
328 const auto ti = AdvectionType::calculateTransmissibilities(problem,
329 element,
330 fvGeometry,
331 curElemVolVars,
332 scvf,
333 elemFluxVarsCache[scvf]);
334
335 // get the rows of the jacobian matrix for the inside/outside scv
336 auto& dI_dJ_inside = A[insideScv.dofIndex()];
337 auto& dI_dJ_outside = A[outsideScv.dofIndex()];
338
339 // precalculate values
340 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
341 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
342 const auto rho_mu_flux_w = rhow_muw*flux_w;
343 const auto rho_mu_flux_n = rhon_mun*flux_n;
344
345 // add the partial derivatives w.r.t all scvs in the element
346 for (const auto& scvJ : scvs(fvGeometry))
347 {
348 const auto globalJ = scvJ.dofIndex();
349 const auto localJ = scvJ.indexInElement();
350
351 // the transmissibily associated with the scvJ
352 const auto tj = ti[localJ];
353
354 // partial derivative of the wetting phase flux w.r.t. p_w
355 const auto tj_up_w = tj*up_w;
356 dI_dJ_inside[globalJ][conti0EqIdx+0][pressureIdx] += tj_up_w;
357 dI_dJ_outside[globalJ][conti0EqIdx+0][pressureIdx] -= tj_up_w;
358
359 // partial derivative of the nonwetting phase flux w.r.t. p_w
360 const auto tj_up_n = tj*up_n;
361 dI_dJ_inside[globalJ][conti0EqIdx+1][pressureIdx] += tj_up_n;
362 dI_dJ_outside[globalJ][conti0EqIdx+1][pressureIdx] -= tj_up_n;
363
364 // partial derivatives w.r.t. S_n (are the negative of those w.r.t sw)
365 // relative permeability contributions only for inside/outside
366 if (localJ == insideScvIdx)
367 {
368 // partial derivative of the wetting phase flux w.r.t. S_n
369 const auto insideSw = insideVolVars.saturation(0);
370 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
371 const auto dFluxW_dSnJ = rho_mu_flux_w*dKrw_dSn_inside*insideWeight_w;
372 dI_dJ_inside[globalJ][conti0EqIdx+0][saturationIdx] += dFluxW_dSnJ;
373 dI_dJ_outside[globalJ][conti0EqIdx+0][saturationIdx] -= dFluxW_dSnJ;
374
375 // partial derivative of the nonwetting phase flux w.r.t. S_n (k_rn contribution)
376 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
377 const auto dFluxN_dSnJ_krn = rho_mu_flux_n*dKrn_dSn_inside*insideWeight_n;
378 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_krn;
379 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_krn;
380
381 // partial derivative of the nonwetting phase flux w.r.t. S_n (p_c contribution)
382 const auto dFluxN_dSnJ_pc = -1.0*tj_up_n*insidefluidMatrixInteraction.dpc_dsw(insideSw);
383 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
384 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
385 }
386 else if (localJ == outsideScvIdx)
387 {
388 // see comments for (globalJ == insideScvIdx)
389 const auto outsideSw = outsideVolVars.saturation(0);
390 const auto dKrw_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrw_dsw(outsideSw);
391 const auto dFluxW_dSnJ = rho_mu_flux_w*dKrw_dSn_outside*outsideWeight_w;
392 dI_dJ_inside[globalJ][conti0EqIdx+0][saturationIdx] += dFluxW_dSnJ;
393 dI_dJ_outside[globalJ][conti0EqIdx+0][saturationIdx] -= dFluxW_dSnJ;
394
395 const auto dKrn_dSn_outside = -1.0*outsidefluidMatrixInteraction.dkrn_dsw(outsideSw);
396 const auto dFluxN_dSnJ_krn = rho_mu_flux_n*dKrn_dSn_outside*outsideWeight_n;
397 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_krn;
398 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_krn;
399
400 const auto dFluxN_dSnJ_pc = -1.0*tj_up_n*outsidefluidMatrixInteraction.dpc_dsw(outsideSw);
401 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
402 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
403 }
404 else
405 {
406 const auto& fluidMatrixInteractionJ
407 = problem.spatialParams().fluidMatrixInteraction(element, scvJ, elemSol);
408
409 const auto swJ = curElemVolVars[scvJ].saturation(0);
410 const auto dFluxN_dSnJ_pc = tj_up_n*fluidMatrixInteractionJ.dpc_dsw(swJ);
411 dI_dJ_inside[globalJ][conti0EqIdx+1][saturationIdx] -= dFluxN_dSnJ_pc;
412 dI_dJ_outside[globalJ][conti0EqIdx+1][saturationIdx] += dFluxN_dSnJ_pc;
413 }
414 }
415 }
416
431 template<class PartialDerivativeMatrices>
432 void addCCDirichletFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
433 const Problem& problem,
434 const Element& element,
435 const FVElementGeometry& fvGeometry,
436 const ElementVolumeVariables& curElemVolVars,
437 const ElementFluxVariablesCache& elemFluxVarsCache,
438 const SubControlVolumeFace& scvf) const
439 {
441
442 // evaluate the current wetting phase Darcy flux and resulting upwind weights
443 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
444 const auto flux_w = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 0, elemFluxVarsCache);
445 const auto flux_n = AdvectionType::flux(problem, element, fvGeometry, curElemVolVars, scvf, 1, elemFluxVarsCache);
446 const auto insideWeight_w = std::signbit(flux_w) ? (1.0 - upwindWeight) : upwindWeight;
447 const auto outsideWeight_w = 1.0 - insideWeight_w;
448 const auto insideWeight_n = std::signbit(flux_n) ? (1.0 - upwindWeight) : upwindWeight;
449 const auto outsideWeight_n = 1.0 - insideWeight_n;
450
451 // get references to the two participating vol vars & parameters
452 const auto insideScvIdx = scvf.insideScvIdx();
453 const auto& insideScv = fvGeometry.scv(insideScvIdx);
454 const auto& insideVolVars = curElemVolVars[insideScvIdx];
455 const auto& outsideVolVars = curElemVolVars[scvf.outsideScvIdx()];
456
457 const auto insidefluidMatrixInteraction = problem.spatialParams().fluidMatrixInteraction(
458 element, insideScv, elementSolution<FVElementGeometry>(insideVolVars.priVars())
459 );
460
461 // some quantities to be reused (rho & mu are constant and thus equal for all cells)
462 static const auto rho_w = insideVolVars.density(0);
463 static const auto rho_n = insideVolVars.density(1);
464 static const auto rhow_muw = rho_w/insideVolVars.viscosity(0);
465 static const auto rhon_mun = rho_n/insideVolVars.viscosity(1);
466 const auto rhowKrw_muw_inside = rho_w*insideVolVars.mobility(0);
467 const auto rhonKrn_mun_inside = rho_n*insideVolVars.mobility(1);
468 const auto rhowKrw_muw_outside = rho_w*outsideVolVars.mobility(0);
469 const auto rhonKrn_mun_outside = rho_n*outsideVolVars.mobility(1);
470
471 // get reference to the inside derivative matrix
472 auto& dI_dI = derivativeMatrices[insideScvIdx];
473
474 // derivative w.r.t. to Sn is the negative of the one w.r.t. Sw
475 const auto insideSw = insideVolVars.saturation(0);
476 const auto dKrw_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrw_dsw(insideSw);
477 const auto dKrn_dSn_inside = -1.0*insidefluidMatrixInteraction.dkrn_dsw(insideSw);
478 const auto dpc_dSn_inside = -1.0*insidefluidMatrixInteraction.dpc_dsw(insideSw);
479
480 const auto tij = elemFluxVarsCache[scvf].advectionTij();
481 // partial derivative of the wetting phase flux w.r.t. p_w
482 const auto up_w = rhowKrw_muw_inside*insideWeight_w + rhowKrw_muw_outside*outsideWeight_w;
483 dI_dI[conti0EqIdx+0][pressureIdx] += tij*up_w;
484
485 // partial derivative of the wetting phase flux w.r.t. S_n
486 dI_dI[conti0EqIdx+0][saturationIdx] += rhow_muw*flux_w*dKrw_dSn_inside*insideWeight_w;
487
488 // partial derivative of the nonwetting phase flux w.r.t. p_w
489 const auto up_n = rhonKrn_mun_inside*insideWeight_n + rhonKrn_mun_outside*outsideWeight_n;
490 dI_dI[conti0EqIdx+1][pressureIdx] += tij*up_n;
491
492 // partial derivative of the nonwetting phase flux w.r.t. S_n (relative permeability derivative contribution)
493 dI_dI[conti0EqIdx+1][saturationIdx] += rhon_mun*flux_n*dKrn_dSn_inside*insideWeight_n;
494
495 // partial derivative of the nonwetting phase flux w.r.t. S_n (capillary pressure derivative contribution)
496 dI_dI[conti0EqIdx+1][saturationIdx] += tij*dpc_dSn_inside*up_n;
497 }
498
510 template<class PartialDerivativeMatrices>
511 void addRobinFluxDerivatives(PartialDerivativeMatrices& derivativeMatrices,
512 const Problem& problem,
513 const Element& element,
514 const FVElementGeometry& fvGeometry,
515 const ElementVolumeVariables& curElemVolVars,
516 const ElementFluxVariablesCache& elemFluxVarsCache,
517 const SubControlVolumeFace& scvf) const
518 { /* TODO maybe forward to problem for the user to implement the Robin derivatives?*/ }
519};
520
521} // end namespace Dumux
522
523#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Element solution classes and factory functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Defines an enumeration for the formulations accepted by the two-phase model.
@ p0s1
first phase pressure and second phase saturation as primary variables
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
auto elementSolution(const Element &element, const SolutionVector &sol, const GridGeometry &gg) -> std::enable_if_t< GridGeometry::discMethod==DiscretizationMethods::box, BoxElementSolution< typename GridGeometry::LocalView, std::decay_t< decltype(std::declval< SolutionVector >()[0])> > >
Make an element solution for box schemes.
Definition: box/elementsolution.hh:118
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr CCTpfa cctpfa
Definition: method.hh:134
constexpr Box box
Definition: method.hh:136
Element-wise calculation of the residual and its derivatives for a two-phase, incompressible test pro...
Definition: 2p/incompressiblelocalresidual.hh:50
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:151
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:92
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:511
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:127
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:270
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:432
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...
Definition: porousmediumflow/immiscible/localresidual.hh:40
Declares all properties used in Dumux.
Element-wise calculation of the residual for problems using the n-phase immiscible fully implicit mod...