2013-05-22 08:16:24 -05:00
|
|
|
/*
|
|
|
|
Copyright 2013 SINTEF ICT, Applied Mathematics.
|
|
|
|
|
|
|
|
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 3 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/>.
|
|
|
|
*/
|
|
|
|
|
|
|
|
#ifndef OPM_GEOPROPS_HEADER_INCLUDED
|
|
|
|
#define OPM_GEOPROPS_HEADER_INCLUDED
|
|
|
|
|
|
|
|
#include <opm/core/grid.h>
|
2014-02-20 06:15:02 -06:00
|
|
|
#include <opm/autodiff/GridHelpers.hpp>
|
2014-07-29 04:40:52 -05:00
|
|
|
#include <opm/core/utility/ErrorMacros.hpp>
|
2014-02-20 06:17:18 -06:00
|
|
|
//#include <opm/core/pressure/tpfa/trans_tpfa.h>
|
|
|
|
#include <opm/core/pressure/tpfa/TransTpfa.hpp>
|
2014-04-15 13:46:45 -05:00
|
|
|
|
2014-07-10 04:43:46 -05:00
|
|
|
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
2015-01-29 02:25:17 -06:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
|
2014-09-20 02:28:25 -05:00
|
|
|
#include <opm/core/utility/platform_dependent/disable_warnings.h>
|
2014-04-15 13:46:45 -05:00
|
|
|
|
2013-05-22 08:16:24 -05:00
|
|
|
#include <Eigen/Eigen>
|
|
|
|
|
2014-07-10 04:53:50 -05:00
|
|
|
#ifdef HAVE_DUNE_CORNERPOINT
|
|
|
|
#include <dune/common/version.hh>
|
|
|
|
#include <dune/grid/CpGrid.hpp>
|
|
|
|
#include <dune/grid/common/mcmgmapper.hh>
|
|
|
|
#endif
|
|
|
|
|
2014-09-20 02:28:25 -05:00
|
|
|
#include <opm/core/utility/platform_dependent/reenable_warnings.h>
|
2014-08-27 12:56:13 -05:00
|
|
|
|
2014-06-27 07:40:30 -05:00
|
|
|
#include <cstddef>
|
2014-04-15 13:46:45 -05:00
|
|
|
|
2013-05-22 08:16:24 -05:00
|
|
|
namespace Opm
|
|
|
|
{
|
|
|
|
|
2013-09-23 06:02:56 -05:00
|
|
|
/// Class containing static geological properties that are
|
2013-09-20 07:55:24 -05:00
|
|
|
/// derived from grid and petrophysical properties:
|
|
|
|
/// - pore volume
|
|
|
|
/// - transmissibilities
|
|
|
|
/// - gravity potentials
|
2013-05-22 08:16:24 -05:00
|
|
|
class DerivedGeology
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
typedef Eigen::ArrayXd Vector;
|
|
|
|
|
|
|
|
/// Construct contained derived geological properties
|
|
|
|
/// from grid and property information.
|
2014-02-20 06:17:18 -06:00
|
|
|
template <class Props, class Grid>
|
2014-07-10 04:43:46 -05:00
|
|
|
DerivedGeology(const Grid& grid,
|
|
|
|
const Props& props ,
|
|
|
|
Opm::EclipseStateConstPtr eclState,
|
2014-12-04 05:29:40 -06:00
|
|
|
const bool use_local_perm,
|
|
|
|
const double* grav = 0
|
|
|
|
|
|
|
|
)
|
2014-02-20 06:15:02 -06:00
|
|
|
: pvol_ (Opm::AutoDiffGrid::numCells(grid))
|
|
|
|
, trans_(Opm::AutoDiffGrid::numFaces(grid))
|
|
|
|
, gpot_ (Vector::Zero(Opm::AutoDiffGrid::cell2Faces(grid).noEntries(), 1))
|
|
|
|
, z_(Opm::AutoDiffGrid::numCells(grid))
|
2013-05-22 08:16:24 -05:00
|
|
|
{
|
2014-07-10 04:53:50 -05:00
|
|
|
int numCells = AutoDiffGrid::numCells(grid);
|
|
|
|
int numFaces = AutoDiffGrid::numFaces(grid);
|
|
|
|
const int *cartDims = AutoDiffGrid::cartDims(grid);
|
|
|
|
int numCartesianCells =
|
|
|
|
cartDims[0]
|
|
|
|
* cartDims[1]
|
|
|
|
* cartDims[2];
|
|
|
|
|
|
|
|
// get the pore volume multipliers from the EclipseState
|
|
|
|
std::vector<double> multpv(numCartesianCells, 1.0);
|
|
|
|
if (eclState->hasDoubleGridProperty("MULTPV")) {
|
|
|
|
multpv = eclState->getDoubleGridProperty("MULTPV")->getData();
|
|
|
|
}
|
|
|
|
|
|
|
|
// get the net-to-gross cell thickness from the EclipseState
|
|
|
|
std::vector<double> ntg(numCartesianCells, 1.0);
|
|
|
|
if (eclState->hasDoubleGridProperty("NTG")) {
|
|
|
|
ntg = eclState->getDoubleGridProperty("NTG")->getData();
|
|
|
|
}
|
|
|
|
|
2015-01-29 02:25:17 -06:00
|
|
|
// get grid from parser.
|
|
|
|
|
|
|
|
// Get original grid cell volume.
|
|
|
|
EclipseGridConstPtr eclgrid = eclState->getEclipseGrid();
|
2015-02-08 21:27:10 -06:00
|
|
|
// Pore volume.
|
|
|
|
// New keywords MINPVF will add some PV due to OPM cpgrid process algorithm.
|
|
|
|
// But the default behavior is to get the comparable pore volume with ECLIPSE.
|
2014-07-10 04:53:50 -05:00
|
|
|
for (int cellIdx = 0; cellIdx < numCells; ++cellIdx) {
|
|
|
|
int cartesianCellIdx = AutoDiffGrid::globalCell(grid)[cellIdx];
|
|
|
|
pvol_[cellIdx] =
|
|
|
|
props.porosity()[cellIdx]
|
|
|
|
* multpv[cartesianCellIdx]
|
2015-02-08 21:27:10 -06:00
|
|
|
* ntg[cartesianCellIdx];
|
2015-02-12 02:16:07 -06:00
|
|
|
if (eclgrid->getMinpvMode() == MinpvMode::ModeEnum::OpmFIL) {
|
2015-02-08 21:27:10 -06:00
|
|
|
pvol_[cellIdx] *= AutoDiffGrid::cellVolume(grid, cellIdx);
|
|
|
|
} else {
|
|
|
|
pvol_[cellIdx] *= eclgrid->getCellVolume(cartesianCellIdx);
|
|
|
|
}
|
2014-07-10 04:53:50 -05:00
|
|
|
}
|
2013-05-22 08:16:24 -05:00
|
|
|
|
|
|
|
// Transmissibility
|
2014-12-01 03:40:50 -06:00
|
|
|
|
2014-07-10 04:53:50 -05:00
|
|
|
Vector htrans(AutoDiffGrid::numCellFaces(grid));
|
2014-02-20 06:17:18 -06:00
|
|
|
Grid* ug = const_cast<Grid*>(& grid);
|
2014-12-04 05:29:40 -06:00
|
|
|
|
2014-12-10 00:34:56 -06:00
|
|
|
if (! use_local_perm) {
|
2014-12-04 05:29:40 -06:00
|
|
|
tpfa_htrans_compute(ug, props.permeability(), htrans.data());
|
2014-12-10 00:34:56 -06:00
|
|
|
}
|
|
|
|
else {
|
2014-12-04 05:29:40 -06:00
|
|
|
tpfa_loc_trans_compute_(grid,props.permeability(),htrans);
|
2014-12-10 00:34:56 -06:00
|
|
|
}
|
2014-07-10 04:53:50 -05:00
|
|
|
|
|
|
|
std::vector<double> mult;
|
|
|
|
multiplyHalfIntersections_(grid, eclState, ntg, htrans, mult);
|
|
|
|
|
|
|
|
// combine the half-face transmissibilites into the final face
|
|
|
|
// transmissibilites.
|
|
|
|
tpfa_trans_compute(ug, htrans.data(), trans_.data());
|
|
|
|
|
|
|
|
// multiply the face transmissibilities with their appropriate
|
|
|
|
// transmissibility multipliers
|
|
|
|
for (int faceIdx = 0; faceIdx < numFaces; faceIdx++) {
|
|
|
|
trans_[faceIdx] *= mult[faceIdx];
|
|
|
|
}
|
2013-05-22 08:16:24 -05:00
|
|
|
|
2013-09-19 07:47:12 -05:00
|
|
|
// Compute z coordinates
|
2014-07-10 04:53:50 -05:00
|
|
|
for (int c = 0; c<numCells; ++c){
|
2015-01-13 13:15:57 -06:00
|
|
|
z_[c] = Opm::UgGridHelpers::cellCentroidCoordinate(grid, c, 2);
|
2013-09-19 07:47:12 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
2013-05-23 07:48:34 -05:00
|
|
|
// Gravity potential
|
|
|
|
std::fill(gravity_, gravity_ + 3, 0.0);
|
2013-05-22 08:16:24 -05:00
|
|
|
if (grav != 0) {
|
2014-07-10 04:53:50 -05:00
|
|
|
const typename Vector::Index nd = AutoDiffGrid::dimensions(grid);
|
|
|
|
typedef typename AutoDiffGrid::ADCell2FacesTraits<Grid>::Type Cell2Faces;
|
|
|
|
Cell2Faces c2f=AutoDiffGrid::cell2Faces(grid);
|
2013-05-22 08:16:24 -05:00
|
|
|
|
2014-06-27 07:40:30 -05:00
|
|
|
std::size_t i = 0;
|
2014-07-10 04:53:50 -05:00
|
|
|
for (typename Vector::Index c = 0; c < numCells; ++c) {
|
|
|
|
const double* const cc = AutoDiffGrid::cellCentroid(grid, c);
|
2013-05-22 08:16:24 -05:00
|
|
|
|
2014-02-20 06:17:18 -06:00
|
|
|
typename Cell2Faces::row_type faces=c2f[c];
|
|
|
|
typedef typename Cell2Faces::row_type::iterator Iter;
|
2013-05-22 08:16:24 -05:00
|
|
|
|
2014-06-27 07:40:30 -05:00
|
|
|
for (Iter f=faces.begin(), end=faces.end(); f!=end; ++f, ++i) {
|
2014-07-10 04:53:50 -05:00
|
|
|
const double* const fc = AutoDiffGrid::faceCentroid(grid, *f);
|
2013-05-22 08:16:24 -05:00
|
|
|
|
|
|
|
for (typename Vector::Index d = 0; d < nd; ++d) {
|
2014-06-27 07:40:30 -05:00
|
|
|
gpot_[i] += grav[d] * (fc[d] - cc[d]);
|
2013-05-22 08:16:24 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
2013-05-23 07:48:34 -05:00
|
|
|
std::copy(grav, grav + nd, gravity_);
|
2013-05-22 08:16:24 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2013-09-24 02:33:33 -05:00
|
|
|
const Vector& poreVolume() const { return pvol_ ;}
|
|
|
|
const Vector& transmissibility() const { return trans_ ;}
|
|
|
|
const Vector& gravityPotential() const { return gpot_ ;}
|
|
|
|
const Vector& z() const { return z_ ;}
|
|
|
|
const double* gravity() const { return gravity_;}
|
2013-05-22 08:16:24 -05:00
|
|
|
|
|
|
|
private:
|
2014-07-10 04:53:50 -05:00
|
|
|
template <class Grid>
|
|
|
|
void multiplyHalfIntersections_(const Grid &grid,
|
|
|
|
Opm::EclipseStateConstPtr eclState,
|
|
|
|
const std::vector<double> &ntg,
|
|
|
|
Vector &halfIntersectTransmissibility,
|
|
|
|
std::vector<double> &intersectionTransMult);
|
|
|
|
|
2014-12-01 03:40:50 -06:00
|
|
|
template <class Grid>
|
|
|
|
void tpfa_loc_trans_compute_(const Grid &grid,
|
|
|
|
const double* perm,
|
|
|
|
Vector &hTrans);
|
|
|
|
|
2013-05-22 08:16:24 -05:00
|
|
|
Vector pvol_ ;
|
|
|
|
Vector trans_;
|
|
|
|
Vector gpot_ ;
|
2013-09-19 07:47:12 -05:00
|
|
|
Vector z_;
|
2013-05-23 07:48:34 -05:00
|
|
|
double gravity_[3]; // Size 3 even if grid is 2-dim.
|
2014-12-01 03:40:50 -06:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2013-05-22 08:16:24 -05:00
|
|
|
};
|
|
|
|
|
2014-12-01 03:40:50 -06:00
|
|
|
|
2015-01-13 13:17:30 -06:00
|
|
|
template <class GridType>
|
|
|
|
inline void DerivedGeology::multiplyHalfIntersections_(const GridType &grid,
|
|
|
|
Opm::EclipseStateConstPtr eclState,
|
|
|
|
const std::vector<double> &ntg,
|
|
|
|
Vector &halfIntersectTransmissibility,
|
|
|
|
std::vector<double> &intersectionTransMult)
|
2014-07-10 04:53:50 -05:00
|
|
|
{
|
2015-01-13 13:17:30 -06:00
|
|
|
int numCells = Opm::AutoDiffGrid::numCells(grid);
|
2014-07-10 04:53:50 -05:00
|
|
|
|
2015-01-13 13:17:30 -06:00
|
|
|
int numIntersections = Opm::AutoDiffGrid::numFaces(grid);
|
2014-07-10 04:53:50 -05:00
|
|
|
intersectionTransMult.resize(numIntersections);
|
|
|
|
std::fill(intersectionTransMult.begin(), intersectionTransMult.end(), 1.0);
|
|
|
|
|
|
|
|
std::shared_ptr<const Opm::TransMult> multipliers = eclState->getTransMult();
|
2015-01-13 13:17:30 -06:00
|
|
|
auto cell2Faces = Opm::UgGridHelpers::cell2Faces(grid);
|
|
|
|
auto faceCells = Opm::AutoDiffGrid::faceCells(grid);
|
|
|
|
const int* global_cell = Opm::UgGridHelpers::globalCell(grid);
|
|
|
|
int cellFaceIdx = 0;
|
2014-07-10 04:53:50 -05:00
|
|
|
|
|
|
|
for (int cellIdx = 0; cellIdx < numCells; ++cellIdx) {
|
|
|
|
// loop over all logically-Cartesian faces of the current cell
|
2015-01-13 13:17:30 -06:00
|
|
|
auto cellFacesRange = cell2Faces[cellIdx];
|
|
|
|
|
|
|
|
for(auto cellFaceIter = cellFacesRange.begin(), cellFaceEnd = cellFacesRange.end();
|
|
|
|
cellFaceIter != cellFaceEnd; ++cellFaceIter, ++cellFaceIdx)
|
2014-07-10 04:53:50 -05:00
|
|
|
{
|
|
|
|
// the index of the current cell in arrays for the logically-Cartesian grid
|
2015-01-13 13:17:30 -06:00
|
|
|
int cartesianCellIdx = global_cell[cellIdx];
|
2014-07-10 04:53:50 -05:00
|
|
|
|
|
|
|
// The index of the face in the compressed grid
|
2015-01-13 13:17:30 -06:00
|
|
|
int faceIdx = *cellFaceIter;
|
2014-07-10 04:53:50 -05:00
|
|
|
|
|
|
|
// the logically-Cartesian direction of the face
|
2015-01-13 13:17:30 -06:00
|
|
|
int faceTag = Opm::UgGridHelpers::faceTag(grid, cellFaceIter);
|
2014-07-10 04:53:50 -05:00
|
|
|
|
|
|
|
// Translate the C face tag into the enum used by opm-parser's TransMult class
|
|
|
|
Opm::FaceDir::DirEnum faceDirection;
|
|
|
|
if (faceTag == 0) // left
|
|
|
|
faceDirection = Opm::FaceDir::XMinus;
|
|
|
|
else if (faceTag == 1) // right
|
|
|
|
faceDirection = Opm::FaceDir::XPlus;
|
|
|
|
else if (faceTag == 2) // back
|
|
|
|
faceDirection = Opm::FaceDir::YMinus;
|
|
|
|
else if (faceTag == 3) // front
|
|
|
|
faceDirection = Opm::FaceDir::YPlus;
|
|
|
|
else if (faceTag == 4) // bottom
|
|
|
|
faceDirection = Opm::FaceDir::ZMinus;
|
|
|
|
else if (faceTag == 5) // top
|
|
|
|
faceDirection = Opm::FaceDir::ZPlus;
|
|
|
|
else
|
|
|
|
OPM_THROW(std::logic_error, "Unhandled face direction: " << faceTag);
|
|
|
|
|
|
|
|
// Account for NTG in horizontal one-sided transmissibilities
|
|
|
|
switch (faceDirection) {
|
|
|
|
case Opm::FaceDir::XMinus:
|
|
|
|
case Opm::FaceDir::XPlus:
|
|
|
|
case Opm::FaceDir::YMinus:
|
|
|
|
case Opm::FaceDir::YPlus:
|
|
|
|
halfIntersectTransmissibility[cellFaceIdx] *= ntg[cartesianCellIdx];
|
|
|
|
break;
|
|
|
|
default:
|
|
|
|
// do nothing for the top and bottom faces
|
|
|
|
break;
|
|
|
|
}
|
|
|
|
|
2014-12-10 00:25:07 -06:00
|
|
|
// Multiplier contribution on this face for MULT[XYZ] logical cartesian multipliers
|
2014-07-10 04:53:50 -05:00
|
|
|
intersectionTransMult[faceIdx] *=
|
|
|
|
multipliers->getMultiplier(cartesianCellIdx, faceDirection);
|
2014-12-10 00:25:07 -06:00
|
|
|
|
|
|
|
// Multiplier contribution on this fase for region multipliers
|
2015-01-13 13:17:30 -06:00
|
|
|
const int cellIdxInside = faceCells(faceIdx, 0);
|
|
|
|
const int cellIdxOutside = faceCells(faceIdx, 1);
|
2014-12-10 00:25:07 -06:00
|
|
|
|
|
|
|
// Do not apply region multipliers in the case of boundary connections
|
|
|
|
if (cellIdxInside < 0 || cellIdxOutside < 0) {
|
|
|
|
continue;
|
|
|
|
}
|
2015-01-13 13:17:30 -06:00
|
|
|
const int cartesianCellIdxInside = global_cell[cellIdxInside];
|
|
|
|
const int cartesianCellIdxOutside = global_cell[cellIdxOutside];
|
2014-12-10 00:25:07 -06:00
|
|
|
// Only apply the region multipliers from the inside
|
|
|
|
if (cartesianCellIdx == cartesianCellIdxInside) {
|
|
|
|
intersectionTransMult[faceIdx] *= multipliers->getRegionMultiplier(cartesianCellIdxInside,cartesianCellIdxOutside,faceDirection);
|
|
|
|
}
|
|
|
|
|
|
|
|
|
2014-07-10 04:53:50 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2014-12-17 03:47:16 -06:00
|
|
|
template <class GridType>
|
2015-01-13 13:17:30 -06:00
|
|
|
inline void DerivedGeology::tpfa_loc_trans_compute_(const GridType& grid,
|
|
|
|
const double* perm,
|
|
|
|
Vector& hTrans){
|
2014-12-01 03:40:50 -06:00
|
|
|
|
|
|
|
// Using Local coordinate system for the transmissibility calculations
|
|
|
|
// hTrans(cellFaceIdx) = K(cellNo,j) * sum( C(:,i) .* N(:,j), 2) / sum(C.*C, 2)
|
|
|
|
// where K is a diagonal permeability tensor, C is the distance from cell centroid
|
|
|
|
// to face centroid and N is the normal vector pointing outwards with norm equal to the face area.
|
|
|
|
// Off-diagonal permeability values are ignored without warning
|
|
|
|
int numCells = AutoDiffGrid::numCells(grid);
|
2015-01-13 13:17:30 -06:00
|
|
|
int cellFaceIdx = 0;
|
|
|
|
auto cell2Faces = Opm::UgGridHelpers::cell2Faces(grid);
|
|
|
|
auto faceCells = Opm::UgGridHelpers::faceCells(grid);
|
|
|
|
|
2014-12-01 03:40:50 -06:00
|
|
|
for (int cellIdx = 0; cellIdx < numCells; ++cellIdx) {
|
|
|
|
// loop over all logically-Cartesian faces of the current cell
|
2015-01-13 13:17:30 -06:00
|
|
|
auto cellFacesRange = cell2Faces[cellIdx];
|
|
|
|
|
|
|
|
for(auto cellFaceIter = cellFacesRange.begin(), cellFaceEnd = cellFacesRange.end();
|
|
|
|
cellFaceIter != cellFaceEnd; ++cellFaceIter, ++cellFaceIdx)
|
2014-12-01 03:40:50 -06:00
|
|
|
{
|
|
|
|
// The index of the face in the compressed grid
|
2015-01-13 13:17:30 -06:00
|
|
|
const int faceIdx = *cellFaceIter;
|
2014-12-01 03:40:50 -06:00
|
|
|
|
|
|
|
// the logically-Cartesian direction of the face
|
2015-01-13 13:17:30 -06:00
|
|
|
const int faceTag = Opm::UgGridHelpers::faceTag(grid, cellFaceIter);
|
2014-12-01 03:40:50 -06:00
|
|
|
|
|
|
|
// d = 0: XPERM d = 4: YPERM d = 8: ZPERM ignores off-diagonal permeability values.
|
|
|
|
const int d = std::floor(faceTag/2) * 4;
|
|
|
|
|
|
|
|
// compute the half transmissibility
|
|
|
|
double dist = 0.0;
|
|
|
|
double cn = 0.0;
|
2015-01-13 13:17:30 -06:00
|
|
|
double sgn = 2.0 * (faceCells(faceIdx, 0) == cellIdx) - 1;
|
|
|
|
const int dim = Opm::UgGridHelpers::dimensions(grid);
|
2014-12-01 03:40:50 -06:00
|
|
|
for (int indx = 0; indx < dim; ++indx) {
|
2015-01-13 13:17:30 -06:00
|
|
|
const double Ci = Opm::UgGridHelpers::faceCentroid(grid, faceIdx)[indx] -
|
|
|
|
Opm::UgGridHelpers::cellCentroidCoordinate(grid, cellIdx, indx);
|
2014-12-01 03:40:50 -06:00
|
|
|
dist += Ci*Ci;
|
2015-01-13 13:17:30 -06:00
|
|
|
cn += sgn * Ci * Opm::UgGridHelpers::faceNormal(grid, faceIdx)[indx];
|
2014-12-01 03:40:50 -06:00
|
|
|
}
|
|
|
|
|
|
|
|
if (cn < 0){
|
|
|
|
switch (d) {
|
|
|
|
case 0:
|
|
|
|
OPM_MESSAGE("Warning: negative X-transmissibility value in cell: " << cellIdx << " replace by absolute value") ;
|
|
|
|
break;
|
|
|
|
case 4:
|
|
|
|
OPM_MESSAGE("Warning: negative Y-transmissibility value in cell: " << cellIdx << " replace by absolute value") ;
|
|
|
|
break;
|
|
|
|
case 8:
|
|
|
|
OPM_MESSAGE("Warning: negative Z-transmissibility value in cell: " << cellIdx << " replace by absolute value") ;
|
|
|
|
break;
|
2014-12-04 05:29:40 -06:00
|
|
|
default:
|
2014-12-10 00:34:56 -06:00
|
|
|
OPM_THROW(std::logic_error, "Inconsistency in the faceTag in cell: " << cellIdx);
|
2014-12-01 03:40:50 -06:00
|
|
|
|
|
|
|
}
|
|
|
|
cn = -cn;
|
|
|
|
}
|
|
|
|
hTrans[cellFaceIdx] = perm[cellIdx*dim*dim + d] * cn / dist;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
}
|
|
|
|
|
2014-07-10 04:53:50 -05:00
|
|
|
}
|
2013-05-22 08:16:24 -05:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
#endif // OPM_GEOPROPS_HEADER_INCLUDED
|