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  Copyright (C) 2008-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
28 #ifndef EWOMS_DARCY_FLUX_MODULE_HH
29 #define EWOMS_DARCY_FLUX_MODULE_HH
30 
33 
34 #include <dune/common/fvector.hh>
35 #include <dune/common/fmatrix.hh>
36 
37 #include <cmath>
38 
39 namespace Ewoms {
40 namespace Properties {
41 NEW_PROP_TAG(MaterialLaw);
42 }
43 
44 template <class TypeTag>
46 
47 template <class TypeTag>
49 
50 template <class TypeTag>
52 
57 template <class TypeTag>
59 {
63 
67  static void registerParameters()
68  { }
69 };
70 
76 template <class TypeTag>
77 class DarcyBaseProblem
78 { };
79 
84 template <class TypeTag>
85 class DarcyIntensiveQuantities
86 {
87  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
88 protected:
89  void update_(const ElementContext &elemCtx, int dofIdx, int timeIdx)
90  { }
91 };
92 
124 template <class TypeTag>
125 class DarcyExtensiveQuantities
126 {
127  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
128  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
129  typedef typename GET_PROP_TYPE(TypeTag, Evaluation) Evaluation;
130  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
131  typedef typename GET_PROP_TYPE(TypeTag, ExtensiveQuantities) Implementation;
132  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
133  typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
134 
135  enum { dimWorld = GridView::dimensionworld };
136  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
137 
138  typedef typename Opm::MathToolbox<Evaluation> Toolbox;
139  typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
140  typedef Dune::FieldVector<Scalar, dimWorld> DimVector;
141  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
142 
143 public:
148  const DimMatrix &intrinsicPermability() const
149  { return K_; }
150 
157  const EvalDimVector &potentialGrad(int phaseIdx) const
158  { return potentialGrad_[phaseIdx]; }
159 
166  const EvalDimVector &filterVelocity(int phaseIdx) const
167  { return filterVelocity_[phaseIdx]; }
168 
178  const Evaluation& volumeFlux(int phaseIdx) const
179  { return volumeFlux_[phaseIdx]; }
180 
181 protected:
182  short upstreamIndex_(int phaseIdx) const
183  { return upstreamDofIdx_[phaseIdx]; }
184 
185  short downstreamIndex_(int phaseIdx) const
186  { return downstreamDofIdx_[phaseIdx]; }
187 
193  void calculateGradients_(const ElementContext &elemCtx,
194  int faceIdx,
195  int timeIdx)
196  {
197  const auto& gradCalc = elemCtx.gradientCalculator();
198  Ewoms::PressureCallback<TypeTag> pressureCallback(elemCtx);
199 
200  const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(faceIdx);
201  const auto &faceNormal = scvf.normal();
202 
203  interiorDofIdx_ = scvf.interiorIndex();
204  exteriorDofIdx_ = scvf.exteriorIndex();
205 
206  // calculate the "raw" pressure gradient
207  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
208  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
209  Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
210  continue;
211  }
212 
213  pressureCallback.setPhaseIndex(phaseIdx);
214  gradCalc.calculateGradient(potentialGrad_[phaseIdx],
215  elemCtx,
216  faceIdx,
217  pressureCallback);
218  Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
219  }
220 
221  // correct the pressure gradients by the gravitational acceleration
222  if (EWOMS_GET_PARAM(TypeTag, bool, EnableGravity)) {
223  // estimate the gravitational acceleration at a given SCV face
224  // using the arithmetic mean
225  const auto& gIn = elemCtx.problem().gravity(elemCtx, interiorDofIdx_, timeIdx);
226  const auto& gEx = elemCtx.problem().gravity(elemCtx, exteriorDofIdx_, timeIdx);
227 
228  const auto &intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
229  const auto &intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
230 
231  const auto &posIn = elemCtx.pos(interiorDofIdx_, timeIdx);
232  const auto &posEx = elemCtx.pos(exteriorDofIdx_, timeIdx);
233  const auto &posFace = scvf.integrationPos();
234 
235  // the distance between the centers of the control volumes
236  DimVector distVecIn(posIn);
237  DimVector distVecEx(posEx);
238  DimVector distVecTotal(posEx);
239 
240  distVecIn -= posFace;
241  distVecEx -= posFace;
242  distVecTotal -= posIn;
243  Scalar absDistTotalSquared = distVecTotal.two_norm2();
244  for (int phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
245  if (!elemCtx.model().phaseIsConsidered(phaseIdx))
246  continue;
247 
248  // calculate the hydrostatic pressure at the integration point of the face
249  auto rhoIn = intQuantsIn.fluidState().density(phaseIdx);
250  auto pStatIn = - rhoIn*(gIn*distVecIn);
251 
252  // the quantities on the exterior side of the face do not influence the
253  // result for the TPFA scheme, so they can be treated as scalar values.
254  Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
255  Scalar pStatEx = - rhoEx*(gEx*distVecEx);
256 
257  // compute the hydrostatic gradient between the two control volumes (this
258  // gradient exhibitis the same direction as the vector between the two
259  // control volume centers and the length (pStaticExterior -
260  // pStaticInterior)/distanceInteriorToExterior
261  EvalDimVector f(distVecTotal);
262  f *= (pStatEx - pStatIn)/absDistTotalSquared;
263 
264  // calculate the final potential gradient
265  potentialGrad_[phaseIdx] += f;
266 
267  for (unsigned i = 0; i < potentialGrad_[phaseIdx].size(); ++i) {
268  if (!std::isfinite(Toolbox::value(potentialGrad_[phaseIdx][i]))) {
269  OPM_THROW(Opm::NumericalProblem,
270  "Non-finite potential gradient for phase '"
271  << FluidSystem::phaseName(phaseIdx) << "'");
272  }
273  }
274  }
275  }
276 
277  Valgrind::SetUndefined(K_);
278  elemCtx.problem().intersectionIntrinsicPermeability(K_, elemCtx, faceIdx, timeIdx);
279  Valgrind::CheckDefined(K_);
280 
281  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
282  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
283  Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
284  continue;
285  }
286 
287  // determine the upstream and downstream DOFs
288  Evaluation tmp = Toolbox::createConstant(0.0);
289  for (unsigned i = 0; i < faceNormal.size(); ++i)
290  tmp += potentialGrad_[phaseIdx][i]*faceNormal[i];
291 
292  if (tmp > 0) {
293  upstreamDofIdx_[phaseIdx] = exteriorDofIdx_;
295  }
296  else {
297  upstreamDofIdx_[phaseIdx] = interiorDofIdx_;
299  }
300 
301  const auto &up = elemCtx.intensiveQuantities(upstreamDofIdx_[phaseIdx], timeIdx);
302  // this is also slightly hacky because it assumes that the derivative of the
303  // flux between two DOFs only depends on the primary variables in the
304  // upstream direction. For non-TPFA flux approximation schemes, this is not
305  // true...
306  if (upstreamDofIdx_[phaseIdx] == interiorDofIdx_)
307  mobility_[phaseIdx] = up.mobility(phaseIdx);
308  else
309  mobility_[phaseIdx] = Toolbox::value(up.mobility(phaseIdx));
310  }
311  }
312 
319  template <class FluidState>
320  void calculateBoundaryGradients_(const ElementContext &elemCtx,
321  int boundaryFaceIdx,
322  int timeIdx,
323  const FluidState& fluidState,
324  const typename FluidSystem::ParameterCache &paramCache)
325  {
326  const auto& gradCalc = elemCtx.gradientCalculator();
327  Ewoms::BoundaryPressureCallback<TypeTag, FluidState> pressureCallback(elemCtx, fluidState);
328 
329  // calculate the pressure gradient
330  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
331  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
332  Valgrind::SetUndefined(potentialGrad_[phaseIdx]);
333  continue;
334  }
335 
336  pressureCallback.setPhaseIndex(phaseIdx);
337  gradCalc.calculateBoundaryGradient(potentialGrad_[phaseIdx],
338  elemCtx,
339  boundaryFaceIdx,
340  pressureCallback);
341  Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
342  }
343 
344  const auto& scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
345  interiorDofIdx_ = scvf.interiorIndex();
346  exteriorDofIdx_ = -1;
347 
348  // calculate the intrinsic permeability
349  const auto &intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
350  K_ = intQuantsIn.intrinsicPermeability();
351 
352  // correct the pressure gradients by the gravitational acceleration
353  if (EWOMS_GET_PARAM(TypeTag, bool, EnableGravity)) {
354  // estimate the gravitational acceleration at a given SCV face
355  // using the arithmetic mean
356  const auto& gIn = elemCtx.problem().gravity(elemCtx, interiorDofIdx_, timeIdx);
357  const auto& posIn = elemCtx.pos(interiorDofIdx_, timeIdx);
358  const auto& posFace = scvf.integrationPos();
359 
360  // the distance between the face center and the center of the control volume
361  DimVector distVecIn(posIn);
362  distVecIn -= posFace;
363  Scalar absDist = distVecIn.two_norm();
364  Scalar gTimesDist = gIn*distVecIn;
365 
366  for (int phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
367  if (!elemCtx.model().phaseIsConsidered(phaseIdx))
368  continue;
369 
370  // calculate the hydrostatic pressure at the integration point of the face
371  Evaluation rhoIn = intQuantsIn.fluidState().density(phaseIdx);
372  Evaluation pStatIn = - gTimesDist*rhoIn;
373 
374  Valgrind::CheckDefined(pStatIn);
375 
376  // compute the hydrostatic gradient between the two control volumes (this
377  // gradient exhibitis the same direction as the vector between the two
378  // control volume centers and the length (pStaticExterior -
379  // pStaticInterior)/distanceInteriorToExterior. Note that for the
380  // boundary, 'pStaticExterior' is zero as the boundary pressure is
381  // defined on boundary face's integration point...
382  EvalDimVector f(distVecIn);
383  f *= - pStatIn/absDist;
384 
385  // calculate the final potential gradient
386  potentialGrad_[phaseIdx] += f;
387 
388  Valgrind::CheckDefined(potentialGrad_[phaseIdx]);
389  for (unsigned i = 0; i < potentialGrad_[phaseIdx].size(); ++i) {
390  if (!std::isfinite(Toolbox::value(potentialGrad_[phaseIdx][i]))) {
391  OPM_THROW(Opm::NumericalProblem,
392  "Non finite potential gradient for phase '"
393  << FluidSystem::phaseName(phaseIdx) << "'");
394  }
395  }
396  }
397  }
398 
399  // determine the upstream and downstream DOFs
400  const auto &faceNormal = scvf.normal();
401 
402  const auto &matParams = elemCtx.problem().materialLawParams(elemCtx, interiorDofIdx_, timeIdx);
403 
404  Scalar kr[numPhases];
405  MaterialLaw::relativePermeabilities(kr, matParams, fluidState);
406  Valgrind::CheckDefined(kr);
407 
408  for (int phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
409  if (!elemCtx.model().phaseIsConsidered(phaseIdx))
410  continue;
411 
412  Evaluation tmp = Toolbox::createConstant(0.0);
413  for (unsigned i = 0; i < faceNormal.size(); ++i)
414  tmp += potentialGrad_[phaseIdx][i]*faceNormal[i];
415 
416  if (tmp > 0) {
417  upstreamDofIdx_[phaseIdx] = exteriorDofIdx_;
419  }
420  else {
421  upstreamDofIdx_[phaseIdx] = interiorDofIdx_;
423  }
424 
425  // take the phase mobility from the DOF in upstream direction
426  if (upstreamDofIdx_[phaseIdx] < 0)
427  mobility_[phaseIdx] =
428  kr[phaseIdx] / FluidSystem::viscosity(fluidState, paramCache, phaseIdx);
429  else
430  mobility_[phaseIdx] = intQuantsIn.mobility(phaseIdx);
431  Valgrind::CheckDefined(mobility_[phaseIdx]);
432  }
433  }
434 
441  void calculateFluxes_(const ElementContext& elemCtx, int scvfIdx, int timeIdx)
442  {
443  const auto &scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
444  const DimVector &normal = scvf.normal();
445  Valgrind::CheckDefined(normal);
446 
447  for (int phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
448  filterVelocity_[phaseIdx] = Toolbox::createConstant(0.0);
449  volumeFlux_[phaseIdx] = Toolbox::createConstant(0.0);
450  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
451  continue;
452  }
453 
454  asImp_().calculateFilterVelocity_(phaseIdx);
455  Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
456  volumeFlux_[phaseIdx] = Toolbox::createConstant(0.0);
457  for (unsigned i = 0; i < normal.size(); ++i)
458  volumeFlux_[phaseIdx] += (filterVelocity_[phaseIdx][i] * normal[i]);
459  }
460  }
461 
468  void calculateBoundaryFluxes_(const ElementContext& elemCtx,
469  int boundaryFaceIdx,
470  int timeIdx)
471  {
472  const auto &scvf = elemCtx.stencil(timeIdx).boundaryFace(boundaryFaceIdx);
473  const DimVector &normal = scvf.normal();
474  Valgrind::CheckDefined(normal);
475 
476  for (int phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
477  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
478  filterVelocity_[phaseIdx] = Toolbox::createConstant(0);
479  volumeFlux_[phaseIdx] = Toolbox::createConstant(0);
480  continue;
481  }
482 
483  asImp_().calculateFilterVelocity_(phaseIdx);
484  Valgrind::CheckDefined(filterVelocity_[phaseIdx]);
485  volumeFlux_[phaseIdx] = Toolbox::createConstant(0.0);
486  for (unsigned i = 0; i < normal.size(); ++i)
487  volumeFlux_[phaseIdx] += (filterVelocity_[phaseIdx][i] * normal[i]);
488  }
489  }
490 
491  void calculateFilterVelocity_(int phaseIdx)
492  {
493 #ifndef NDEBUG
494  assert(std::isfinite(Toolbox::value(mobility_[phaseIdx])));
495  for (unsigned i = 0; i < K_.M(); ++ i)
496  for (unsigned j = 0; j < K_.N(); ++ j)
497  assert(std::isfinite(K_[i][j]));
498 #endif
499 
500  K_.mv(potentialGrad_[phaseIdx], filterVelocity_[phaseIdx]);
501  filterVelocity_[phaseIdx] *= - mobility_[phaseIdx];
502 
503 #ifndef NDEBUG
504  for (unsigned i = 0; i < filterVelocity_[phaseIdx].size(); ++ i)
505  assert(std::isfinite(Toolbox::value(filterVelocity_[phaseIdx][i])));
506 #endif
507  }
508 
509 private:
510  Implementation &asImp_()
511  { return *static_cast<Implementation*>(this); }
512 
513  const Implementation &asImp_() const
514  { return *static_cast<const Implementation*>(this); }
515 
516 protected:
517  // intrinsic permeability tensor and its square root
518  DimMatrix K_;
519 
520  // interior, exterior, upstream and downstream DOFs
523  short upstreamDofIdx_[numPhases];
524  short downstreamDofIdx_[numPhases];
525 
526  // mobilities of all fluid phases [1 / (Pa s)]
527  Evaluation mobility_[numPhases];
528 
529  // filter velocities of all phases [m/s]
530  EvalDimVector filterVelocity_[numPhases];
531 
532  // the volumetric flux of all fluid phases over the control
533  // volume's face [m^3/s / m^2]
534  Evaluation volumeFlux_[numPhases];
535 
536  // pressure potential gradients of all phases [Pa / m]
537  EvalDimVector potentialGrad_[numPhases];
538 };
539 
540 } // namespace Ewoms
541 
542 #endif
const EvalDimVector & filterVelocity(int phaseIdx) const
Return the filter velocity of a fluid phase at the face's integration point [m/s].
Definition: darcyfluxmodule.hh:166
short interiorDofIdx_
Definition: darcyfluxmodule.hh:521
void setPhaseIndex(short phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:155
const DimMatrix & intrinsicPermability() const
Returns the intrinsic permeability tensor for a given sub-control volume face.
Definition: darcyfluxmodule.hh:148
Specifies a flux module which uses the Darcy relation.
Definition: darcyfluxmodule.hh:58
void calculateFluxes_(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
Calculate the volumetric fluxes of all phases.
Definition: darcyfluxmodule.hh:441
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
short upstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:523
const EvalDimVector & potentialGrad(int phaseIdx) const
Return the pressure potential gradient of a fluid phase at the face's integration point [Pa/m]...
Definition: darcyfluxmodule.hh:157
void update_(const ElementContext &elemCtx, int dofIdx, int timeIdx)
Definition: darcyfluxmodule.hh:89
const Evaluation & volumeFlux(int phaseIdx) const
Return the volume flux of a fluid phase at the face's integration point .
Definition: darcyfluxmodule.hh:178
DarcyIntensiveQuantities< TypeTag > FluxIntensiveQuantities
Definition: darcyfluxmodule.hh:60
NEW_PROP_TAG(Grid)
The type of the DUNE grid.
Provides the defaults for the parameters required by the Darcy velocity approach. ...
Definition: darcyfluxmodule.hh:51
void setPhaseIndex(short phaseIdx)
Set the index of the fluid phase for which the pressure should be returned.
Definition: quantitycallbacks.hh:101
short upstreamIndex_(int phaseIdx) const
Definition: darcyfluxmodule.hh:182
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:125
short downstreamIndex_(int phaseIdx) const
Definition: darcyfluxmodule.hh:185
Evaluation mobility_[numPhases]
Definition: darcyfluxmodule.hh:527
void calculateBoundaryGradients_(const ElementContext &elemCtx, int boundaryFaceIdx, int timeIdx, const FluidState &fluidState, const typename FluidSystem::ParameterCache &paramCache)
Calculate the gradients at the grid boundary which are required to determine the volumetric fluxes...
Definition: darcyfluxmodule.hh:320
Definition: baseauxiliarymodule.hh:35
Callback class for a phase pressure.
Definition: quantitycallbacks.hh:78
void calculateBoundaryFluxes_(const ElementContext &elemCtx, int boundaryFaceIdx, int timeIdx)
Calculate the volumetric fluxes at a boundary face of all fluid phases.
Definition: darcyfluxmodule.hh:468
Evaluation volumeFlux_[numPhases]
Definition: darcyfluxmodule.hh:534
DimMatrix K_
Definition: darcyfluxmodule.hh:518
Provides the Darcy flux module.
Definition: darcyfluxmodule.hh:48
static void registerParameters()
Register all run-time parameters for the flux module.
Definition: darcyfluxmodule.hh:67
void calculateFilterVelocity_(int phaseIdx)
Definition: darcyfluxmodule.hh:491
This method contains all callback classes for quantities that are required by some extensive quantiti...
Provides the intensive quantities for the Darcy flux module.
Definition: darcyfluxmodule.hh:45
EvalDimVector filterVelocity_[numPhases]
Definition: darcyfluxmodule.hh:530
DarcyExtensiveQuantities< TypeTag > FluxExtensiveQuantities
Definition: darcyfluxmodule.hh:61
void calculateGradients_(const ElementContext &elemCtx, int faceIdx, int timeIdx)
Calculate the gradients which are required to determine the volumetric fluxes.
Definition: darcyfluxmodule.hh:193
short downstreamDofIdx_[numPhases]
Definition: darcyfluxmodule.hh:524
Defines the common properties required by the porous medium multi-phase models.
EvalDimVector potentialGrad_[numPhases]
Definition: darcyfluxmodule.hh:537
DarcyBaseProblem< TypeTag > FluxBaseProblem
Definition: darcyfluxmodule.hh:62
#define EWOMS_GET_PARAM(TypeTag, ParamType, ParamName)
Retrieve a runtime parameter.
Definition: parametersystem.hh:95
short exteriorDofIdx_
Definition: darcyfluxmodule.hh:522