darcyfluxmodule.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 This file is part of the Open Porous Media project (OPM).
5
6 OPM 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 2 of the License, or
9 (at your option) any later version.
10
11 OPM 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 OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
30#ifndef EWOMS_DARCY_FLUX_MODULE_HH
31#define EWOMS_DARCY_FLUX_MODULE_HH
32
33#include <dune/common/fvector.hh>
34#include <dune/common/fmatrix.hh>
35
36#include <opm/common/Exceptions.hpp>
37
38#include <opm/material/common/Valgrind.hpp>
39
43
44namespace Opm {
45
46template <class TypeTag>
47class DarcyIntensiveQuantities;
48
49template <class TypeTag>
50class DarcyExtensiveQuantities;
51
52template <class TypeTag>
53class DarcyBaseProblem;
54
59template <class TypeTag>
61{
65
69 static void registerParameters()
70 { }
71};
72
78template <class TypeTag>
80{ };
81
86template <class TypeTag>
88{
90protected:
91 void update_(const ElementContext&,
92 unsigned,
93 unsigned)
94 { }
95};
96
113template <class TypeTag>
115{
123
124 enum { dimWorld = GridView::dimensionworld };
125 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
126
127 using Toolbox = MathToolbox<Evaluation>;
128 using ParameterCache = typename FluidSystem::template ParameterCache<Evaluation>;
129 using EvalDimVector = Dune::FieldVector<Evaluation, dimWorld>;
130 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
131 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
132
133public:
138 const DimMatrix& intrinsicPermability() const
139 { return K_; }
140
147 const EvalDimVector& potentialGrad(unsigned phaseIdx) const
148 { return potentialGrad_[phaseIdx]; }
149
156 const EvalDimVector& filterVelocity(unsigned phaseIdx) const
157 { return filterVelocity_[phaseIdx]; }
158
168 const Evaluation& volumeFlux(unsigned phaseIdx) const
169 { return volumeFlux_[phaseIdx]; }
170
171protected:
172 short upstreamIndex_(unsigned phaseIdx) const
173 { return upstreamDofIdx_[phaseIdx]; }
174
175 short downstreamIndex_(unsigned phaseIdx) const
176 { return downstreamDofIdx_[phaseIdx]; }
177
183 void calculateGradients_(const ElementContext& elemCtx,
184 unsigned faceIdx,
185 unsigned timeIdx)
186 {
187 const auto& gradCalc = elemCtx.gradientCalculator();
188 PressureCallback<TypeTag> pressureCallback(elemCtx);
189
190 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
191 const auto& faceNormal = scvf.normal();
192
193 unsigned i = scvf.interiorIndex();
194 unsigned j = scvf.exteriorIndex();
195 interiorDofIdx_ = static_cast<short>(i);
196 exteriorDofIdx_ = static_cast<short>(j);
197 unsigned focusDofIdx = elemCtx.focusDofIndex();
198
199 // calculate the "raw" pressure gradient
200 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
201 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
202 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
203 continue;
204 }
205
206 pressureCallback.setPhaseIndex(phaseIdx);
207 gradCalc.calculateGradient(potentialGrad_[phaseIdx],
208 elemCtx,
209 faceIdx,
210 pressureCallback);
211 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
212 }
213
214 // correct the pressure gradients by the gravitational acceleration
215 if (Parameters::Get<Parameters::EnableGravity>()) {
216 // estimate the gravitational acceleration at a given SCV face
217 // using the arithmetic mean
218 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
219 const auto& gEx = elemCtx.problem().gravity(elemCtx, j, timeIdx);
220
221 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
222 const auto& intQuantsEx = elemCtx.intensiveQuantities(j, timeIdx);
223
224 const auto& posIn = elemCtx.pos(i, timeIdx);
225 const auto& posEx = elemCtx.pos(j, timeIdx);
226 const auto& posFace = scvf.integrationPos();
227
228 // the distance between the centers of the control volumes
229 DimVector distVecIn(posIn);
230 DimVector distVecEx(posEx);
231 DimVector distVecTotal(posEx);
232
233 distVecIn -= posFace;
234 distVecEx -= posFace;
235 distVecTotal -= posIn;
236 Scalar absDistTotalSquared = distVecTotal.two_norm2();
237 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
238 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
239 continue;
240
241 // calculate the hydrostatic pressure at the integration point of the face
242 Evaluation pStatIn;
243
244 if (std::is_same<Scalar, Evaluation>::value ||
245 interiorDofIdx_ == static_cast<int>(focusDofIdx))
246 {
247 const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
248 pStatIn = - rhoIn*(gIn*distVecIn);
249 }
250 else {
251 Scalar rhoIn = Toolbox::value(intQuantsIn.fluidState().density(phaseIdx));
252 pStatIn = - rhoIn*(gIn*distVecIn);
253 }
254
255 // the quantities on the exterior side of the face do not influence the
256 // result for the TPFA scheme, so they can be treated as scalar values.
257 Evaluation pStatEx;
258
259 if (std::is_same<Scalar, Evaluation>::value ||
260 exteriorDofIdx_ == static_cast<int>(focusDofIdx))
261 {
262 const Evaluation& rhoEx = intQuantsEx.fluidState().density(phaseIdx);
263 pStatEx = - rhoEx*(gEx*distVecEx);
264 }
265 else {
266 Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
267 pStatEx = - rhoEx*(gEx*distVecEx);
268 }
269
270 // compute the hydrostatic gradient between the two control volumes (this
271 // gradient exhibitis the same direction as the vector between the two
272 // control volume centers and the length (pStaticExterior -
273 // pStaticInterior)/distanceInteriorToExterior
274 Dune::FieldVector<Evaluation, dimWorld> f(distVecTotal);
275 f *= (pStatEx - pStatIn)/absDistTotalSquared;
276
277 // calculate the final potential gradient
278 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
279 potentialGrad_[phaseIdx][dimIdx] += f[dimIdx];
280
281 for (unsigned dimIdx = 0; dimIdx < potentialGrad_[phaseIdx].size(); ++dimIdx) {
282 if (!isfinite(potentialGrad_[phaseIdx][dimIdx])) {
283 throw NumericalProblem("Non-finite potential gradient for phase '"
284 + std::string(FluidSystem::phaseName(phaseIdx))+"'");
285 }
286 }
287 }
288 }
289
290 Valgrind::SetUndefined(K_);
291 elemCtx.problem().intersectionIntrinsicPermeability(K_, elemCtx, faceIdx, timeIdx);
292 Valgrind::CheckDefined(K_);
293
294 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
295 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
296 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
297 continue;
298 }
299
300 // determine the upstream and downstream DOFs
301 Evaluation tmp = 0.0;
302 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
303 tmp += potentialGrad_[phaseIdx][dimIdx]*faceNormal[dimIdx];
304
305 if (tmp > 0) {
308 }
309 else {
312 }
313
314 // we only carry the derivatives along if the upstream DOF is the one which
315 // we currently focus on
316 const auto& up = elemCtx.intensiveQuantities(upstreamDofIdx_[phaseIdx], timeIdx);
317 if (upstreamDofIdx_[phaseIdx] == static_cast<int>(focusDofIdx))
318 mobility_[phaseIdx] = up.mobility(phaseIdx);
319 else
320 mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
321 }
322 }
323
330 template <class FluidState>
331 void calculateBoundaryGradients_(const ElementContext& elemCtx,
332 unsigned boundaryFaceIdx,
333 unsigned timeIdx,
334 const FluidState& fluidState)
335 {
336 const auto& gradCalc = elemCtx.gradientCalculator();
337 BoundaryPressureCallback<TypeTag, FluidState> pressureCallback(elemCtx, fluidState);
338
339 // calculate the pressure gradient
340 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
341 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
342 Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
343 continue;
344 }
345
346 pressureCallback.setPhaseIndex(phaseIdx);
347 gradCalc.calculateBoundaryGradient(potentialGrad_[phaseIdx],
348 elemCtx,
349 boundaryFaceIdx,
350 pressureCallback);
351 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
352 }
353
354 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
355 auto i = scvf.interiorIndex();
356 interiorDofIdx_ = static_cast<short>(i);
357 exteriorDofIdx_ = -1;
358 int focusDofIdx = elemCtx.focusDofIndex();
359
360 // calculate the intrinsic permeability
361 const auto& intQuantsIn = elemCtx.intensiveQuantities(i, timeIdx);
362 K_ = intQuantsIn.intrinsicPermeability();
363
364 // correct the pressure gradients by the gravitational acceleration
365 if (Parameters::Get<Parameters::EnableGravity>()) {
366 // estimate the gravitational acceleration at a given SCV face
367 // using the arithmetic mean
368 const auto& gIn = elemCtx.problem().gravity(elemCtx, i, timeIdx);
369 const auto& posIn = elemCtx.pos(i, timeIdx);
370 const auto& posFace = scvf.integrationPos();
371
372 // the distance between the face center and the center of the control volume
373 DimVector distVecIn(posIn);
374 distVecIn -= posFace;
375 Scalar absDistSquared = distVecIn.two_norm2();
376 Scalar gTimesDist = gIn*distVecIn;
377
378 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
379 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
380 continue;
381
382 // calculate the hydrostatic pressure at the integration point of the face
383 Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
384 Evaluation pStatIn = - gTimesDist*rhoIn;
385
386 Valgrind::CheckDefined(pStatIn);
387 // compute the hydrostatic gradient between the control volume and face integration
388 // point. This gradient exhibitis the same direction as the vector between the
389 // control volume center and face integration point (-distVecIn / absDist) and the
390 // length of the vector is -pStaticIn / absDist. Note that the two negatives become
391 // + and the 1 / (absDist * absDist) -> absDistSquared.
392 EvalDimVector f(distVecIn);
393 f *= pStatIn / absDistSquared;
394
395 // calculate the final potential gradient
396 for (unsigned dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
397 potentialGrad_[phaseIdx][dimIdx] += f[dimIdx];
398
399 Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
400 for (unsigned dimIdx = 0; dimIdx < potentialGrad_[phaseIdx].size(); ++dimIdx) {
401 if (!isfinite(potentialGrad_[phaseIdx][dimIdx])) {
402 throw NumericalProblem("Non-finite potential gradient for phase '"
403 + std::string(FluidSystem::phaseName(phaseIdx))+"'");
404 }
405 }
406 }
407 }
408
409 // determine the upstream and downstream DOFs
410 const auto& faceNormal = scvf.normal();
411
412 const auto& matParams = elemCtx.problem().materialLawParams(elemCtx, i, timeIdx);
413
414 Scalar kr[numPhases];
415 MaterialLaw::relativePermeabilities(kr, matParams, fluidState);
416 Valgrind::CheckDefined(kr);
417
418 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
419 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
420 continue;
421
422 Evaluation tmp = 0.0;
423 for (unsigned dimIdx = 0; dimIdx < faceNormal.size(); ++dimIdx)
424 tmp += potentialGrad_[phaseIdx][dimIdx]*faceNormal[dimIdx];
425
426 if (tmp > 0) {
429 }
430 else {
433 }
434
435 // take the phase mobility from the DOF in upstream direction
436 if (upstreamDofIdx_[phaseIdx] < 0) {
437 if (interiorDofIdx_ == focusDofIdx)
438 mobility_[phaseIdx] =
439 kr[phaseIdx] / fluidState.viscosity(phaseIdx);
440 else
441 mobility_[phaseIdx] =
442 Toolbox::value(kr[phaseIdx])
443 / Toolbox::value(fluidState.viscosity(phaseIdx));
444 }
445 else if (upstreamDofIdx_[phaseIdx] != focusDofIdx)
446 mobility_[phaseIdx] = Toolbox::value(intQuantsIn.mobility(phaseIdx));
447 else
448 mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
449 Valgrind::CheckDefined(mobility_[phaseIdx]);
450 }
451 }
452
459 void calculateFluxes_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
460 {
461 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
462 const DimVector& normal = scvf.normal();
463 Valgrind::CheckDefined(normal);
464
465 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
466 filterVelocity_[phaseIdx] = 0.0;
467 volumeFlux_[phaseIdx] = 0.0;
468 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
469 continue;
470
471 asImp_().calculateFilterVelocity_(phaseIdx);
472 Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
473
474 volumeFlux_[phaseIdx] = 0.0;
475 for (unsigned i = 0; i < normal.size(); ++i)
476 volumeFlux_[phaseIdx] += filterVelocity_[phaseIdx][i] * normal[i];
477 }
478 }
479
486 void calculateBoundaryFluxes_(const ElementContext& elemCtx,
487 unsigned boundaryFaceIdx,
488 unsigned timeIdx)
489 {
490 const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
491 const DimVector& normal = scvf.normal();
492 Valgrind::CheckDefined(normal);
493
494 for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
495 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
496 filterVelocity_[phaseIdx] = 0.0;
497 volumeFlux_[phaseIdx] = 0.0;
498 continue;
499 }
500
501 asImp_().calculateFilterVelocity_(phaseIdx);
502 Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
503 volumeFlux_[phaseIdx] = 0.0;
504 for (unsigned i = 0; i < normal.size(); ++i)
505 volumeFlux_[phaseIdx] += filterVelocity_[phaseIdx][i] * normal[i];
506 }
507 }
508
509 void calculateFilterVelocity_(unsigned phaseIdx)
510 {
511#ifndef NDEBUG
512 assert(isfinite(mobility_[phaseIdx]));
513 for (unsigned i = 0; i < K_.M(); ++ i)
514 for (unsigned j = 0; j < K_.N(); ++ j)
515 assert(std::isfinite(K_[i][j]));
516#endif
517
518 K_.mv(potentialGrad_[phaseIdx], filterVelocity_[phaseIdx]);
519 filterVelocity_[phaseIdx] *= - mobility_[phaseIdx];
520
521#ifndef NDEBUG
522 for (unsigned i = 0; i < filterVelocity_[phaseIdx].size(); ++ i)
523 assert(isfinite(filterVelocity_[phaseIdx][i]));
524#endif
525 }
526
527private:
528 Implementation& asImp_()
529 { return *static_cast<Implementation*>(this); }
530
531 const Implementation& asImp_() const
532 { return *static_cast<const Implementation*>(this); }
533
534protected:
535 // intrinsic permeability tensor and its square root
536 DimMatrix K_;
537
538 // mobilities of all fluid phases [1 / (Pa s)]
539 Evaluation mobility_[numPhases];
540
541 // filter velocities of all phases [m/s]
542 EvalDimVector filterVelocity_[numPhases];
543
544 // the volumetric flux of all fluid phases over the control
545 // volume's face [m^3/s / m^2]
546 Evaluation volumeFlux_[numPhases];
547
548 // pressure potential gradients of all phases [Pa / m]
549 EvalDimVector potentialGrad_[numPhases];
550
551 // upstream, downstream, interior and exterior DOFs
552 short upstreamDofIdx_[numPhases];
553 short downstreamDofIdx_[numPhases];
556};
557
558} // namespace Opm
559
560#endif
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:133
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:163
Provides the defaults for the parameters required by the Darcy velocity approach.
Definition: darcyfluxmodule.hh:80
Provides the Darcy flux module.
Definition: darcyfluxmodule.hh:115
Evaluation mobility_[numPhases]
Definition: darcyfluxmodule.hh:539
void calculateBoundaryGradients_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx, const FluidState &fluidState)
Calculate the gradients at the grid boundary which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:331
short interiorDofIdx_
Definition: darcyfluxmodule.hh:554
const EvalDimVector & filterVelocity(unsigned phaseIdx) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: darcyfluxmodule.hh:156
void calculateBoundaryFluxes_(const ElementContext &elemCtx, unsigned boundaryFaceIdx, unsigned timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition: darcyfluxmodule.hh:486
short upstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:172
EvalDimVector potentialGrad_[numPhases]
Definition: darcyfluxmodule.hh:549
EvalDimVector filterVelocity_[numPhases]
Definition: darcyfluxmodule.hh:542
const EvalDimVector & potentialGrad(unsigned phaseIdx) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m].
Definition: darcyfluxmodule.hh:147
short upstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:552
void calculateFluxes_(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: darcyfluxmodule.hh:459
short downstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:553
const Evaluation & volumeFlux(unsigned phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: darcyfluxmodule.hh:168
short downstreamIndex_(unsigned phaseIdx) const
Definition: darcyfluxmodule.hh:175
Evaluation volumeFlux_[numPhases]
Definition: darcyfluxmodule.hh:546
void calculateFilterVelocity_(unsigned phaseIdx)
Definition: darcyfluxmodule.hh:509
void calculateGradients_(const ElementContext &elemCtx, unsigned faceIdx, unsigned timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:183
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:555
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:138
DimMatrix K_
Definition: darcyfluxmodule.hh:536
Provides the intensive quantities for the Darcy flux module.
Definition: darcyfluxmodule.hh:88
void update_(const ElementContext &, unsigned, unsigned)
Definition: darcyfluxmodule.hh:91
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:84
void setPhaseIndex(unsigned phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:108
Defines the common parameters for the porous medium multi-phase models.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:37
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:235
This method contains all callback classes for quantities that are required by some extensive quantiti...
Specifies a flux module which uses the Darcy relation.
Definition: darcyfluxmodule.hh:61
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: darcyfluxmodule.hh:69