1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511
|
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 2 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
Consult the COPYING file in the top-level source directory of this
module for the precise wording of the license and the list of
copyright holders.
*/
/*!
* \file
*
* \brief This file contains the flux module that uses transmissibilities
*
* The transmissibility approach to fluxes used here is limited
* to the two-point flux approximation
*/
#ifndef EWOMS_TRANS_FLUX_MODULE_HH
#define EWOMS_TRANS_FLUX_MODULE_HH
#include "multiphasebaseproperties.hh"
#include <opm/models/utils/signum.hh>
#include <opm/material/common/Valgrind.hpp>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
namespace Opm {
template <class TypeTag>
class TransIntensiveQuantities;
template <class TypeTag>
class TransExtensiveQuantities;
template <class TypeTag>
class TransBaseProblem;
/*!
* \brief Specifies a flux module which uses transmissibilities.
*/
template <class TypeTag>
struct TransFluxModule
{
using FluxIntensiveQuantities = TransIntensiveQuantities<TypeTag>;
using FluxExtensiveQuantities = TransExtensiveQuantities<TypeTag>;
using FluxBaseProblem = TransBaseProblem<TypeTag>;
/*!
* \brief Register all run-time parameters for the flux module.
*/
static void registerParameters()
{ }
};
/*!
* \brief Provides the defaults for the parameters required by the
* transmissibility based volume flux calculation.
*/
template <class TypeTag>
class TransBaseProblem
{ };
/*!
* \brief Provides the intensive quantities for the transmissibility based flux module
*/
template <class TypeTag>
class TransIntensiveQuantities
{
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
protected:
void update_(const ElementContext&, unsigned, unsigned)
{ }
};
/*!
* \brief Provides the transmissibility based flux module
*/
template <class TypeTag>
class TransExtensiveQuantities
{
using Implementation = GetPropType<TypeTag, Properties::ExtensiveQuantities>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using Scalar = GetPropType<TypeTag, Properties::Scalar>;
using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
using GridView = GetPropType<TypeTag, Properties::GridView>;
using MaterialLaw = GetPropType<TypeTag, Properties::MaterialLaw>;
using Discretization = GetPropType<TypeTag, Properties::Discretization>;
enum { dimWorld = GridView::dimensionworld };
enum { numPhases = FluidSystem::numPhases };
typedef MathToolbox<Evaluation> Toolbox;
typedef Dune::FieldVector<Scalar, dimWorld> DimVector;
typedef Dune::FieldVector<Evaluation, dimWorld> EvalDimVector;
typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
public:
/*!
* \brief Return the intrinsic permeability tensor at a face [m^2]
*/
const DimMatrix& intrinsicPermeability() const
{
throw std::logic_error("The ECL transmissibility module does not provide an explicit intrinsic permeability");
}
/*!
* \brief Return the pressure potential gradient of a fluid phase at the
* face's integration point [Pa/m]
*
* \param phaseIdx The index of the fluid phase
*/
const EvalDimVector& potentialGrad(unsigned) const
{
throw std::logic_error("The ECL transmissibility module does not provide explicit potential gradients");
}
/*!
* \brief Return the gravity corrected pressure difference between the interior and
* the exterior of a face.
*
* \param phaseIdx The index of the fluid phase
*/
const Evaluation& pressureDifference(unsigned phaseIdx) const
{ return pressureDifference_[phaseIdx]; }
/*!
* \brief Return the filter velocity of a fluid phase at the face's integration point
* [m/s]
*
* \param phaseIdx The index of the fluid phase
*/
const EvalDimVector& filterVelocity(unsigned) const
{
throw std::logic_error("The ECL transmissibility module does not provide explicit filter velocities");
}
/*!
* \brief Return the volume flux of a fluid phase at the face's integration point
* \f$[m^3/s / m^2]\f$
*
* This is the fluid volume of a phase per second and per square meter of face
* area.
*
* \param phaseIdx The index of the fluid phase
*/
const Evaluation& volumeFlux(unsigned phaseIdx) const
{ return volumeFlux_[phaseIdx]; }
protected:
/*!
* \brief Returns the local index of the degree of freedom in which is
* in upstream direction.
*
* i.e., the DOF which exhibits a higher effective pressure for
* the given phase.
*/
unsigned upstreamIndex_(unsigned phaseIdx) const
{
assert(phaseIdx < numPhases);
return upIdx_[phaseIdx];
}
/*!
* \brief Returns the local index of the degree of freedom in which is
* in downstream direction.
*
* i.e., the DOF which exhibits a lower effective pressure for the
* given phase.
*/
unsigned downstreamIndex_(unsigned phaseIdx) const
{
assert(phaseIdx < numPhases);
return dnIdx_[phaseIdx];
}
void updateSolvent(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
{ asImp_().updateVolumeFluxTrans(elemCtx, scvfIdx, timeIdx); }
void updatePolymer(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
{ asImp_().updateShearMultipliers(elemCtx, scvfIdx, timeIdx); }
/*!
* \brief Update the required gradients for interior faces
*/
void calculateGradients_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
{
Valgrind::SetUndefined(*this);
// only valied for element center finite volume discretization
static const bool isEcfv = std::is_same<Discretization, EcfvDiscretization<TypeTag> >::value;
static_assert(isEcfv);
const auto& stencil = elemCtx.stencil(timeIdx);
const auto& scvf = stencil.interiorFace(scvfIdx);
interiorDofIdx_ = scvf.interiorIndex();
exteriorDofIdx_ = scvf.exteriorIndex();
assert(interiorDofIdx_ != exteriorDofIdx_);
unsigned I = stencil.globalSpaceIndex(interiorDofIdx_);
unsigned J = stencil.globalSpaceIndex(exteriorDofIdx_);
Scalar trans = transmissibility_(elemCtx, scvfIdx, timeIdx);
// estimate the gravity correction: for performance reasons we use a simplified
// approach for this flux module that assumes that gravity is constant and always
// acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
const auto& intQuantsEx = elemCtx.intensiveQuantities(exteriorDofIdx_, timeIdx);
Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
Scalar zEx = dofCenterDepth_(elemCtx, exteriorDofIdx_, timeIdx);
// the distances from the DOF's depths. (i.e., the additional depth of the
// exterior DOF)
Scalar distZ = zIn - zEx;
for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
if (!FluidSystem::phaseIsActive(phaseIdx))
continue;
// check shortcut: if the mobility of the phase is zero in the interior as
// well as the exterior DOF, we can skip looking at the phase.
if (intQuantsIn.mobility(phaseIdx) <= 0.0 &&
intQuantsEx.mobility(phaseIdx) <= 0.0)
{
upIdx_[phaseIdx] = interiorDofIdx_;
dnIdx_[phaseIdx] = exteriorDofIdx_;
pressureDifference_[phaseIdx] = 0.0;
volumeFlux_[phaseIdx] = 0.0;
continue;
}
// do the gravity correction: compute the hydrostatic pressure for the
// external at the depth of the internal one
const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
Scalar rhoEx = Toolbox::value(intQuantsEx.fluidState().density(phaseIdx));
Evaluation rhoAvg = (rhoIn + rhoEx)/2;
const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
Evaluation pressureExterior = Toolbox::value(intQuantsEx.fluidState().pressure(phaseIdx));
pressureExterior += rhoAvg*(distZ*g);
pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
// decide the upstream index for the phase. for this we make sure that the
// degree of freedom which is regarded upstream if both pressures are equal
// is always the same: if the pressure is equal, the DOF with the lower
// global index is regarded to be the upstream one.
if (pressureDifference_[phaseIdx] > 0.0) {
upIdx_[phaseIdx] = exteriorDofIdx_;
dnIdx_[phaseIdx] = interiorDofIdx_;
}
else if (pressureDifference_[phaseIdx] < 0.0) {
upIdx_[phaseIdx] = interiorDofIdx_;
dnIdx_[phaseIdx] = exteriorDofIdx_;
}
else {
// if the pressure difference is zero, we chose the DOF which has the
// larger volume associated to it as upstream DOF
Scalar Vin = elemCtx.dofVolume(interiorDofIdx_, /*timeIdx=*/0);
Scalar Vex = elemCtx.dofVolume(exteriorDofIdx_, /*timeIdx=*/0);
if (Vin > Vex) {
upIdx_[phaseIdx] = interiorDofIdx_;
dnIdx_[phaseIdx] = exteriorDofIdx_;
}
else if (Vin < Vex) {
upIdx_[phaseIdx] = exteriorDofIdx_;
dnIdx_[phaseIdx] = interiorDofIdx_;
}
else {
assert(Vin == Vex);
// if the volumes are also equal, we pick the DOF which exhibits the
// smaller global index
if (I < J) {
upIdx_[phaseIdx] = interiorDofIdx_;
dnIdx_[phaseIdx] = exteriorDofIdx_;
}
else {
upIdx_[phaseIdx] = exteriorDofIdx_;
dnIdx_[phaseIdx] = interiorDofIdx_;
}
}
}
// this is slightly hacky because in the automatic differentiation case, it
// only works for the element centered finite volume method. for ebos this
// does not matter, though.
unsigned upstreamIdx = upstreamIndex_(phaseIdx);
const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
if (upstreamIdx == interiorDofIdx_)
volumeFlux_[phaseIdx] =
pressureDifference_[phaseIdx]*up.mobility(phaseIdx)*(-trans);
else
volumeFlux_[phaseIdx] =
pressureDifference_[phaseIdx]*(Toolbox::value(up.mobility(phaseIdx))*(-trans));
}
}
/*!
* \brief Update the required gradients for boundary faces
*/
template <class FluidState>
void calculateBoundaryGradients_(const ElementContext& elemCtx,
unsigned scvfIdx,
unsigned timeIdx,
const FluidState& exFluidState)
{
const auto& stencil = elemCtx.stencil(timeIdx);
const auto& scvf = stencil.boundaryFace(scvfIdx);
interiorDofIdx_ = scvf.interiorIndex();
Scalar trans = transmissibilityBoundary_(elemCtx, scvfIdx, timeIdx);
// estimate the gravity correction: for performance reasons we use a simplified
// approach for this flux module that assumes that gravity is constant and always
// acts into the downwards direction. (i.e., no centrifuge experiments, sorry.)
Scalar g = elemCtx.problem().gravity()[dimWorld - 1];
const auto& intQuantsIn = elemCtx.intensiveQuantities(interiorDofIdx_, timeIdx);
// this is quite hacky because the dune grid interface does not provide a
// cellCenterDepth() method (so we ask the problem to provide it). The "good"
// solution would be to take the Z coordinate of the element centroids, but since
// ECL seems to like to be inconsistent on that front, it needs to be done like
// here...
Scalar zIn = dofCenterDepth_(elemCtx, interiorDofIdx_, timeIdx);
Scalar zEx = scvf.integrationPos()[dimWorld - 1];
// the distances from the DOF's depths. (i.e., the additional depth of the
// exterior DOF)
Scalar distZ = zIn - zEx;
for (unsigned phaseIdx=0; phaseIdx < numPhases; phaseIdx++) {
if (!FluidSystem::phaseIsActive(phaseIdx))
continue;
// do the gravity correction: compute the hydrostatic pressure for the
// integration position
const Evaluation& rhoIn = intQuantsIn.fluidState().density(phaseIdx);
const auto& rhoEx = exFluidState.density(phaseIdx);
Evaluation rhoAvg = (rhoIn + rhoEx)/2;
const Evaluation& pressureInterior = intQuantsIn.fluidState().pressure(phaseIdx);
Evaluation pressureExterior = exFluidState.pressure(phaseIdx);
pressureExterior += rhoAvg*(distZ*g);
pressureDifference_[phaseIdx] = pressureExterior - pressureInterior;
// decide the upstream index for the phase. for this we make sure that the
// degree of freedom which is regarded upstream if both pressures are equal
// is always the same: if the pressure is equal, the DOF with the lower
// global index is regarded to be the upstream one.
if (pressureDifference_[phaseIdx] > 0.0) {
upIdx_[phaseIdx] = -1;
dnIdx_[phaseIdx] = interiorDofIdx_;
}
else {
upIdx_[phaseIdx] = interiorDofIdx_;
dnIdx_[phaseIdx] = -1;
}
short upstreamIdx = upstreamIndex_(phaseIdx);
if (upstreamIdx == interiorDofIdx_) {
// this is slightly hacky because in the automatic differentiation case, it
// only works for the element centered finite volume method. for ebos this
// does not matter, though.
const auto& up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
volumeFlux_[phaseIdx] =
pressureDifference_[phaseIdx]*up.mobility(phaseIdx)*(-trans);
}
else {
// compute the phase mobility using the material law parameters of the
// interior element. \todo {this could probably be done more efficiently}
const auto& matParams =
elemCtx.problem().materialLawParams(elemCtx,
interiorDofIdx_,
/*timeIdx=*/0);
std::array<typename FluidState::Scalar,numPhases> kr;
MaterialLaw::relativePermeabilities(kr, matParams, exFluidState);
const auto& mob = kr[phaseIdx]/exFluidState.viscosity(phaseIdx);
volumeFlux_[phaseIdx] =
pressureDifference_[phaseIdx]*mob*(-trans);
}
}
}
/*!
* \brief Update the volumetric fluxes for all fluid phases on the interior faces of the context
*/
void calculateFluxes_(const ElementContext&, unsigned, unsigned)
{ }
void calculateBoundaryFluxes_(const ElementContext&, unsigned, unsigned)
{}
private:
Scalar transmissibility_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx) const
{
const auto& stencil = elemCtx.stencil(timeIdx);
const auto& face = stencil.interiorFace(scvfIdx);
const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
const auto& exteriorPos = stencil.subControlVolume(face.exteriorIndex()).globalPos();
auto distVec0 = face.integrationPos() - interiorPos;
auto distVec1 = face.integrationPos() - exteriorPos;
Scalar ndotDistIn = std::abs(face.normal() * distVec0);
Scalar ndotDistExt = std::abs(face.normal() * distVec1);
Scalar distSquaredIn = distVec0 * distVec0;
Scalar distSquaredExt = distVec1 * distVec1;
const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
const auto& K1mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.exteriorIndex(), timeIdx);
// the permeability per definition aligns with the grid
// we only support diagonal permeability tensor
// and can therefore neglect off-diagonal values
int idx = 0;
Scalar val = 0.0;
for (unsigned i = 0; i < dimWorld; ++ i){
if (std::abs(face.normal()[i]) > val) {
val = std::abs(face.normal()[i]);
idx = i;
}
}
const Scalar& K0 = K0mat[idx][idx];
const Scalar& K1 = K1mat[idx][idx];
const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
const Scalar T1 = K1 * ndotDistExt / distSquaredExt;
return T0 * T1 / (T0 + T1);
}
Scalar transmissibilityBoundary_(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx) const
{
const auto& stencil = elemCtx.stencil(timeIdx);
const auto& face = stencil.interiorFace(scvfIdx);
const auto& interiorPos = stencil.subControlVolume(face.interiorIndex()).globalPos();
auto distVec0 = face.integrationPos() - interiorPos;
Scalar ndotDistIn = face.normal() * distVec0;
Scalar distSquaredIn = distVec0 * distVec0;
const auto& K0mat = elemCtx.problem().intrinsicPermeability(elemCtx, face.interiorIndex(), timeIdx);
// the permeability per definition aligns with the grid
// we only support diagonal permeability tensor
// and can therefore neglect off-diagonal values
int idx = 0;
Scalar val = 0.0;
for (unsigned i = 0; i < dimWorld; ++ i){
if (std::abs(face.normal()[i]) > val) {
val = std::abs(face.normal()[i]);
idx = i;
}
}
const Scalar& K0 = K0mat[idx][idx];
const Scalar T0 = K0 * ndotDistIn / distSquaredIn;
return T0;
}
template <class Context>
Scalar dofCenterDepth_(const Context& context, unsigned spaceIdx, unsigned timeIdx) const
{
const auto& pos = context.pos(spaceIdx, timeIdx);
return pos[dimWorld-1];
}
Implementation& asImp_()
{ return *static_cast<Implementation*>(this); }
const Implementation& asImp_() const
{ return *static_cast<const Implementation*>(this); }
// the volumetric flux of all phases [m^3/s]
Evaluation volumeFlux_[numPhases];
// the difference in effective pressure between the exterior and the interior degree
// of freedom [Pa]
Evaluation pressureDifference_[numPhases];
// the local indices of the interior and exterior degrees of freedom
unsigned short interiorDofIdx_;
unsigned short exteriorDofIdx_;
short upIdx_[numPhases];
short dnIdx_[numPhases];
};
} // namespace Opm
#endif
|