mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Add support for NNC in the simulator
1) NNC are added the grad, div and average operators 2) NNC are added the upwindSelector 3) NNC transmissibilities are added to the face transmissibilities
This commit is contained in:
@@ -1,5 +1,6 @@
|
|||||||
/*
|
/*
|
||||||
Copyright 2013 SINTEF ICT, Applied Mathematics.
|
Copyright 2013 SINTEF ICT, Applied Mathematics.
|
||||||
|
Copyright 2015 IRIS
|
||||||
|
|
||||||
This file is part of the Open Porous Media project (OPM).
|
This file is part of the Open Porous Media project (OPM).
|
||||||
|
|
||||||
@@ -24,7 +25,8 @@
|
|||||||
#include <opm/autodiff/GridHelpers.hpp>
|
#include <opm/autodiff/GridHelpers.hpp>
|
||||||
#include <opm/core/grid.h>
|
#include <opm/core/grid.h>
|
||||||
#include <opm/core/utility/ErrorMacros.hpp>
|
#include <opm/core/utility/ErrorMacros.hpp>
|
||||||
|
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
|
||||||
|
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
||||||
#include <iostream>
|
#include <iostream>
|
||||||
#include <vector>
|
#include <vector>
|
||||||
|
|
||||||
@@ -58,6 +60,109 @@ struct HelperOps
|
|||||||
/// Extract for each cell the sum of all its adjacent faces' (signed) values.
|
/// Extract for each cell the sum of all its adjacent faces' (signed) values.
|
||||||
M fulldiv;
|
M fulldiv;
|
||||||
|
|
||||||
|
/// Non-neighboring connections
|
||||||
|
typedef Eigen::Array<int, Eigen::Dynamic, 2, Eigen::RowMajor> TwoColInt;
|
||||||
|
TwoColInt nnc_cells;
|
||||||
|
|
||||||
|
/// The NNC transmissibilities
|
||||||
|
V nnc_trans;
|
||||||
|
|
||||||
|
/// Constructs all helper vectors and matrices.
|
||||||
|
template<class Grid>
|
||||||
|
HelperOps(const Grid& grid, Opm::EclipseStateConstPtr eclState)
|
||||||
|
{
|
||||||
|
using namespace AutoDiffGrid;
|
||||||
|
const int nc = numCells(grid);
|
||||||
|
const int nf = numFaces(grid);
|
||||||
|
// Define some neighbourhood-derived helper arrays.
|
||||||
|
|
||||||
|
TwoColInt nbi;
|
||||||
|
extractInternalFaces(grid, internal_faces, nbi);
|
||||||
|
int num_internal=internal_faces.size();
|
||||||
|
// num_connections may also include non-neighboring connections
|
||||||
|
int num_connections = num_internal;
|
||||||
|
int numNNC = 0;
|
||||||
|
|
||||||
|
// handle non-neighboring connections
|
||||||
|
std::shared_ptr<const NNC> nnc = eclState->getNNC();
|
||||||
|
if (nnc->hasNNC()) {
|
||||||
|
numNNC = nnc->numNNC();
|
||||||
|
num_connections += numNNC;
|
||||||
|
//std::cout << "Added " << numNNC << " NNC" <<std::endl;
|
||||||
|
nbi.resize(num_internal, 2);
|
||||||
|
|
||||||
|
// the nnc's acts on global indicies and must be mapped to cell indicies
|
||||||
|
size_t cartesianSize = eclState->getEclipseGrid()->getCartesianSize();
|
||||||
|
std::vector<int> global2localIdx(cartesianSize,0);
|
||||||
|
for (int i = 0; i< nc; ++i) {
|
||||||
|
global2localIdx[grid.global_cell[i]] = i;
|
||||||
|
}
|
||||||
|
const std::vector<size_t>& NNC1 = nnc->nnc1();
|
||||||
|
const std::vector<size_t>& NNC2 = nnc->nnc2();
|
||||||
|
nnc_cells.resize(numNNC,2);
|
||||||
|
for (int i = 0; i < numNNC; ++i) {
|
||||||
|
nnc_cells(i,0) = global2localIdx[NNC1[i]];
|
||||||
|
nnc_cells(i,1) = global2localIdx[NNC2[i]];
|
||||||
|
}
|
||||||
|
// store the nnc transmissibilities for later usage.
|
||||||
|
nnc_trans = Eigen::Map<const V>(nnc->trans().data(), numNNC);
|
||||||
|
} else {
|
||||||
|
nnc_trans.resize(0);
|
||||||
|
nnc_cells.resize(0,0);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
// std::cout << "nbi = \n" << nbi << std::endl;
|
||||||
|
// Create matrices.
|
||||||
|
ngrad.resize(num_connections, nc);
|
||||||
|
caver.resize(num_connections, nc);
|
||||||
|
typedef Eigen::Triplet<double> Tri;
|
||||||
|
std::vector<Tri> ngrad_tri;
|
||||||
|
std::vector<Tri> caver_tri;
|
||||||
|
ngrad_tri.reserve(2*num_connections);
|
||||||
|
caver_tri.reserve(2*num_connections);
|
||||||
|
for (int i = 0; i < num_internal; ++i) {
|
||||||
|
ngrad_tri.emplace_back(i, nbi(i,0), 1.0);
|
||||||
|
ngrad_tri.emplace_back(i, nbi(i,1), -1.0);
|
||||||
|
caver_tri.emplace_back(i, nbi(i,0), 0.5);
|
||||||
|
caver_tri.emplace_back(i, nbi(i,1), 0.5);
|
||||||
|
}
|
||||||
|
// add contribution from NNC
|
||||||
|
if (nnc->hasNNC()) {
|
||||||
|
for (int i = 0; i < numNNC; ++i) {
|
||||||
|
ngrad_tri.emplace_back(i+num_internal, nnc_cells(i,0), 1.0);
|
||||||
|
ngrad_tri.emplace_back(i+num_internal, nnc_cells(i,1), -1.0);
|
||||||
|
caver_tri.emplace_back(i+num_internal, nnc_cells(i,0), 0.5);
|
||||||
|
caver_tri.emplace_back(i+num_internal, nnc_cells(i,1), 0.5);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
ngrad.setFromTriplets(ngrad_tri.begin(), ngrad_tri.end());
|
||||||
|
caver.setFromTriplets(caver_tri.begin(), caver_tri.end());
|
||||||
|
grad = -ngrad;
|
||||||
|
div = ngrad.transpose();
|
||||||
|
|
||||||
|
std::vector<Tri> fullngrad_tri;
|
||||||
|
fullngrad_tri.reserve(2*(nf+numNNC));
|
||||||
|
typename ADFaceCellTraits<Grid>::Type nb=faceCellsToEigen(grid);
|
||||||
|
for (int i = 0; i < nf; ++i) {
|
||||||
|
if (nb(i,0) >= 0) {
|
||||||
|
fullngrad_tri.emplace_back(i, nb(i,0), 1.0);
|
||||||
|
}
|
||||||
|
if (nb(i,1) >= 0) {
|
||||||
|
fullngrad_tri.emplace_back(i, nb(i,1), -1.0);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
// add contribution from NNC
|
||||||
|
if (nnc->hasNNC()) {
|
||||||
|
for (int i = 0; i < numNNC; ++i) {
|
||||||
|
fullngrad_tri.emplace_back(i+nf, nnc_cells(i,0), 1.0);
|
||||||
|
fullngrad_tri.emplace_back(i+nf, nnc_cells(i,1), -1.0);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
fullngrad.resize(nf+numNNC, nc);
|
||||||
|
fullngrad.setFromTriplets(fullngrad_tri.begin(), fullngrad_tri.end());
|
||||||
|
fulldiv = fullngrad.transpose();
|
||||||
|
}
|
||||||
/// Constructs all helper vectors and matrices.
|
/// Constructs all helper vectors and matrices.
|
||||||
template<class Grid>
|
template<class Grid>
|
||||||
HelperOps(const Grid& grid)
|
HelperOps(const Grid& grid)
|
||||||
@@ -66,7 +171,7 @@ struct HelperOps
|
|||||||
const int nc = numCells(grid);
|
const int nc = numCells(grid);
|
||||||
const int nf = numFaces(grid);
|
const int nf = numFaces(grid);
|
||||||
// Define some neighbourhood-derived helper arrays.
|
// Define some neighbourhood-derived helper arrays.
|
||||||
typedef Eigen::Array<int, Eigen::Dynamic, 2, Eigen::RowMajor> TwoColInt;
|
//typedef Eigen::Array<int, Eigen::Dynamic, 2, Eigen::RowMajor> TwoColInt;
|
||||||
TwoColInt nbi;
|
TwoColInt nbi;
|
||||||
extractInternalFaces(grid, internal_faces, nbi);
|
extractInternalFaces(grid, internal_faces, nbi);
|
||||||
int num_internal=internal_faces.size();
|
int num_internal=internal_faces.size();
|
||||||
@@ -127,11 +232,15 @@ struct HelperOps
|
|||||||
const IFIndex nif = h.internal_faces.size();
|
const IFIndex nif = h.internal_faces.size();
|
||||||
typename ADFaceCellTraits<Grid>::Type
|
typename ADFaceCellTraits<Grid>::Type
|
||||||
face_cells = faceCellsToEigen(g);
|
face_cells = faceCellsToEigen(g);
|
||||||
assert(nif == ifaceflux.size());
|
|
||||||
|
// num connections may possibly include NNCs
|
||||||
|
size_t num_nnc = h.nnc_trans.size();
|
||||||
|
size_t num_connections= nif + num_nnc;
|
||||||
|
assert(num_connections == ifaceflux.size());
|
||||||
|
|
||||||
// Define selector structure.
|
// Define selector structure.
|
||||||
typedef typename Eigen::Triplet<Scalar> Triplet;
|
typedef typename Eigen::Triplet<Scalar> Triplet;
|
||||||
std::vector<Triplet> s; s.reserve(nif);
|
std::vector<Triplet> s; s.reserve(num_connections);
|
||||||
for (IFIndex iface = 0; iface < nif; ++iface) {
|
for (IFIndex iface = 0; iface < nif; ++iface) {
|
||||||
const int f = h.internal_faces[iface];
|
const int f = h.internal_faces[iface];
|
||||||
const int c1 = face_cells(f,0);
|
const int c1 = face_cells(f,0);
|
||||||
@@ -144,9 +253,15 @@ struct HelperOps
|
|||||||
|
|
||||||
s.push_back(Triplet(iface, c, Scalar(1)));
|
s.push_back(Triplet(iface, c, Scalar(1)));
|
||||||
}
|
}
|
||||||
|
if (num_nnc>0) {
|
||||||
|
for (int i = 0; i<num_nnc; ++i) {
|
||||||
|
const int c = (ifaceflux[i+nif] >= 0) ? h.nnc_cells(i,0) : h.nnc_cells(i,1);
|
||||||
|
s.push_back(Triplet(i+nif,c,Scalar(1)));
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
// Assemble explicit selector operator.
|
// Assemble explicit selector operator.
|
||||||
select_.resize(nif, numCells(g));
|
select_.resize(num_connections, numCells(g));
|
||||||
select_.setFromTriplets(s.begin(), s.end());
|
select_.setFromTriplets(s.begin(), s.end());
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@@ -55,6 +55,7 @@ namespace Opm {
|
|||||||
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
||||||
/// \param[in] wells well structure
|
/// \param[in] wells well structure
|
||||||
/// \param[in] linsolver linear solver
|
/// \param[in] linsolver linear solver
|
||||||
|
/// \param[in] eclState eclipse state
|
||||||
/// \param[in] has_disgas turn on dissolved gas
|
/// \param[in] has_disgas turn on dissolved gas
|
||||||
/// \param[in] has_vapoil turn on vaporized oil feature
|
/// \param[in] has_vapoil turn on vaporized oil feature
|
||||||
/// \param[in] terminal_output request output to cout/cerr
|
/// \param[in] terminal_output request output to cout/cerr
|
||||||
@@ -65,11 +66,12 @@ namespace Opm {
|
|||||||
const RockCompressibility* rock_comp_props,
|
const RockCompressibility* rock_comp_props,
|
||||||
const Wells* wells,
|
const Wells* wells,
|
||||||
const NewtonIterationBlackoilInterface& linsolver,
|
const NewtonIterationBlackoilInterface& linsolver,
|
||||||
|
Opm::EclipseStateConstPtr eclState,
|
||||||
const bool has_disgas,
|
const bool has_disgas,
|
||||||
const bool has_vapoil,
|
const bool has_vapoil,
|
||||||
const bool terminal_output)
|
const bool terminal_output)
|
||||||
: Base(param, grid, fluid, geo, rock_comp_props, wells, linsolver,
|
: Base(param, grid, fluid, geo, rock_comp_props, wells, linsolver,
|
||||||
has_disgas, has_vapoil, terminal_output)
|
eclState, has_disgas, has_vapoil, terminal_output)
|
||||||
{
|
{
|
||||||
}
|
}
|
||||||
};
|
};
|
||||||
|
|||||||
@@ -31,6 +31,7 @@
|
|||||||
#include <opm/autodiff/LinearisedBlackoilResidual.hpp>
|
#include <opm/autodiff/LinearisedBlackoilResidual.hpp>
|
||||||
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
|
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
|
||||||
#include <opm/autodiff/BlackoilModelEnums.hpp>
|
#include <opm/autodiff/BlackoilModelEnums.hpp>
|
||||||
|
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
|
||||||
|
|
||||||
#include <array>
|
#include <array>
|
||||||
|
|
||||||
@@ -117,6 +118,7 @@ namespace Opm {
|
|||||||
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
||||||
/// \param[in] wells well structure
|
/// \param[in] wells well structure
|
||||||
/// \param[in] linsolver linear solver
|
/// \param[in] linsolver linear solver
|
||||||
|
/// \param[in] eclState eclipse state
|
||||||
/// \param[in] has_disgas turn on dissolved gas
|
/// \param[in] has_disgas turn on dissolved gas
|
||||||
/// \param[in] has_vapoil turn on vaporized oil feature
|
/// \param[in] has_vapoil turn on vaporized oil feature
|
||||||
/// \param[in] terminal_output request output to cout/cerr
|
/// \param[in] terminal_output request output to cout/cerr
|
||||||
@@ -127,6 +129,7 @@ namespace Opm {
|
|||||||
const RockCompressibility* rock_comp_props,
|
const RockCompressibility* rock_comp_props,
|
||||||
const Wells* wells,
|
const Wells* wells,
|
||||||
const NewtonIterationBlackoilInterface& linsolver,
|
const NewtonIterationBlackoilInterface& linsolver,
|
||||||
|
Opm::EclipseStateConstPtr eclState,
|
||||||
const bool has_disgas,
|
const bool has_disgas,
|
||||||
const bool has_vapoil,
|
const bool has_vapoil,
|
||||||
const bool terminal_output);
|
const bool terminal_output);
|
||||||
|
|||||||
@@ -42,6 +42,7 @@
|
|||||||
#include <opm/core/utility/Units.hpp>
|
#include <opm/core/utility/Units.hpp>
|
||||||
#include <opm/core/well_controls.h>
|
#include <opm/core/well_controls.h>
|
||||||
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
||||||
|
#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
|
||||||
|
|
||||||
#include <cassert>
|
#include <cassert>
|
||||||
#include <cmath>
|
#include <cmath>
|
||||||
@@ -144,6 +145,7 @@ namespace detail {
|
|||||||
const RockCompressibility* rock_comp_props,
|
const RockCompressibility* rock_comp_props,
|
||||||
const Wells* wells,
|
const Wells* wells,
|
||||||
const NewtonIterationBlackoilInterface& linsolver,
|
const NewtonIterationBlackoilInterface& linsolver,
|
||||||
|
Opm::EclipseStateConstPtr eclState,
|
||||||
const bool has_disgas,
|
const bool has_disgas,
|
||||||
const bool has_vapoil,
|
const bool has_vapoil,
|
||||||
const bool terminal_output)
|
const bool terminal_output)
|
||||||
@@ -156,7 +158,7 @@ namespace detail {
|
|||||||
, active_(detail::activePhases(fluid.phaseUsage()))
|
, active_(detail::activePhases(fluid.phaseUsage()))
|
||||||
, canph_ (detail::active2Canonical(fluid.phaseUsage()))
|
, canph_ (detail::active2Canonical(fluid.phaseUsage()))
|
||||||
, cells_ (detail::buildAllCells(Opm::AutoDiffGrid::numCells(grid)))
|
, cells_ (detail::buildAllCells(Opm::AutoDiffGrid::numCells(grid)))
|
||||||
, ops_ (grid)
|
, ops_ (grid, eclState)
|
||||||
, wops_ (wells_)
|
, wops_ (wells_)
|
||||||
, has_disgas_(has_disgas)
|
, has_disgas_(has_disgas)
|
||||||
, has_vapoil_(has_vapoil)
|
, has_vapoil_(has_vapoil)
|
||||||
@@ -837,9 +839,13 @@ namespace detail {
|
|||||||
// Set up the common parts of the mass balance equations
|
// Set up the common parts of the mass balance equations
|
||||||
// for each active phase.
|
// for each active phase.
|
||||||
const V transi = subset(geo_.transmissibility(), ops_.internal_faces);
|
const V transi = subset(geo_.transmissibility(), ops_.internal_faces);
|
||||||
|
const V trans_nnc = ops_.nnc_trans;
|
||||||
|
V trans_all(transi.size() + trans_nnc.size());
|
||||||
|
trans_all << transi, trans_nnc;
|
||||||
|
|
||||||
const std::vector<ADB> kr = computeRelPerm(state);
|
const std::vector<ADB> kr = computeRelPerm(state);
|
||||||
for (int phaseIdx = 0; phaseIdx < fluid_.numPhases(); ++phaseIdx) {
|
for (int phaseIdx = 0; phaseIdx < fluid_.numPhases(); ++phaseIdx) {
|
||||||
asImpl().computeMassFlux(phaseIdx, transi, kr[canph_[phaseIdx]], state.canonical_phase_pressures[canph_[phaseIdx]], state);
|
asImpl().computeMassFlux(phaseIdx, trans_all, kr[canph_[phaseIdx]], state.canonical_phase_pressures[canph_[phaseIdx]], state);
|
||||||
|
|
||||||
residual_.material_balance_eq[ phaseIdx ] =
|
residual_.material_balance_eq[ phaseIdx ] =
|
||||||
pvdt_ * (rq_[phaseIdx].accum[1] - rq_[phaseIdx].accum[0])
|
pvdt_ * (rq_[phaseIdx].accum[1] - rq_[phaseIdx].accum[0])
|
||||||
|
|||||||
@@ -341,6 +341,7 @@ namespace Opm
|
|||||||
rock_comp_props_,
|
rock_comp_props_,
|
||||||
wells,
|
wells,
|
||||||
solver_,
|
solver_,
|
||||||
|
eclipse_state_,
|
||||||
has_disgas_,
|
has_disgas_,
|
||||||
has_vapoil_,
|
has_vapoil_,
|
||||||
terminal_output_));
|
terminal_output_));
|
||||||
|
|||||||
Reference in New Issue
Block a user