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