Make Connecting Cell's Background Porosity Available in ScheduleGrid

The porosity is needed for WDFACCOR-related connection output to the
restart file.
This commit is contained in:
Bård Skaflestad 2023-11-10 12:52:49 +01:00
parent 73c58dba48
commit f790e81e9a
12 changed files with 215 additions and 123 deletions

View File

@ -16,40 +16,42 @@
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef COMPLETED_CELLS
#define COMPLETED_CELLS
#include <optional>
#include <unordered_map>
#include <opm/input/eclipse/EclipseState/Grid/GridDims.hpp>
#include <array>
#include <cstddef>
#include <optional>
#include <unordered_map>
#include <utility>
namespace Opm {
class CompletedCells {
class CompletedCells
{
public:
struct Cell {
struct Cell
{
std::size_t global_index;
std::size_t i, j, k;
struct Props{
std::size_t active_index;
double permx;
double permy;
double permz;
int satnum;
int pvtnum;
double ntg;
struct Props
{
std::size_t active_index{};
double permx{};
double permy{};
double permz{};
double poro{};
int satnum{};
int pvtnum{};
double ntg{};
bool operator==(const Props& other) const{
return this->active_index == other.active_index &&
this->permx == other.permx &&
this->permy == other.permy &&
this->permz == other.permz &&
this->satnum == other.satnum &&
this->pvtnum == other.pvtnum &&
this->ntg == other.ntg;
}
bool operator==(const Props& other) const;
static Props serializationTestObject();
template<class Serializer>
void serializeOp(Serializer& serializer)
@ -57,46 +59,23 @@ public:
serializer(this->permx);
serializer(this->permy);
serializer(this->permz);
serializer(this->poro);
serializer(this->satnum);
serializer(this->pvtnum);
serializer(this->ntg);
}
static Props serializationTestObject(){
Props props;
props.permx = 10.0;
props.permy = 78.0;
props.permz = 45.4;
props.satnum = 3;
props.pvtnum = 5;
props.ntg = 45.1;
return props;
}
};
std::optional<Props> props;
std::size_t active_index() const;
bool is_active() const;
double depth;
std::array<double, 3> dimensions;
double depth{};
std::array<double, 3> dimensions{};
bool operator==(const Cell& other) const {
return this->global_index == other.global_index &&
this->i == other.i &&
this->j == other.j &&
this->k == other.k &&
this->depth == other.depth &&
this->dimensions == other.dimensions &&
this->props == other.props;
}
bool operator==(const Cell& other) const;
static Cell serializationTestObject() {
Cell cell(0,1,1,1);
cell.depth = 12345;
cell.dimensions = {1.0,2.0,3.0};
return cell;
}
static Cell serializationTestObject();
template<class Serializer>
void serializeOp(Serializer& serializer)
@ -123,6 +102,7 @@ public:
CompletedCells() = default;
explicit CompletedCells(const GridDims& dims);
CompletedCells(std::size_t nx, std::size_t ny, std::size_t nz);
const Cell& get(std::size_t i, std::size_t j, std::size_t k) const;
std::pair<bool, Cell&> try_get(std::size_t i, std::size_t j, std::size_t k);
@ -141,5 +121,5 @@ private:
std::unordered_map<std::size_t, Cell> cells;
};
}
#endif
#endif // COMPLETED_CELLS

View File

@ -21,17 +21,30 @@
#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
#include <cstddef>
#include <string>
namespace Opm {
class EclipseGrid;
class FieldPropsManager;
class ScheduleGrid {
} // namespace Opm
namespace Opm {
class ScheduleGrid
{
public:
ScheduleGrid(const EclipseGrid& ecl_grid, const FieldPropsManager& fpm, CompletedCells& completed_cells);
ScheduleGrid(const EclipseGrid& ecl_grid,
const FieldPropsManager& fpm,
CompletedCells& completed_cells);
explicit ScheduleGrid(CompletedCells& completed_cells);
const CompletedCells::Cell& get_cell(std::size_t i, std::size_t j, std::size_t k) const;
const CompletedCells::Cell&
get_cell(std::size_t i, std::size_t j, std::size_t k) const;
const Opm::EclipseGrid* get_grid() const;
private:
@ -40,8 +53,6 @@ private:
CompletedCells& cells;
};
} // namespace Opm
}
#endif
#endif // SCHEDULE_GRID

View File

@ -19,50 +19,105 @@
#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
bool Opm::CompletedCells::Cell::Props::operator==(const Props& other) const
{
return (this->active_index == other.active_index)
&& (this->permx == other.permx)
&& (this->permy == other.permy)
&& (this->permz == other.permz)
&& (this->poro == other.poro)
&& (this->satnum == other.satnum)
&& (this->pvtnum == other.pvtnum)
&& (this->ntg == other.ntg)
;
}
Opm::CompletedCells::Cell::Props
Opm::CompletedCells::Cell::Props::serializationTestObject()
{
Props props;
props.permx = 10.0;
props.permy = 78.0;
props.permz = 45.4;
props.poro = 0.321;
props.satnum = 3;
props.pvtnum = 5;
props.ntg = 45.1;
return props;
}
bool Opm::CompletedCells::Cell::operator==(const Cell& other) const
{
return (this->global_index == other.global_index)
&& (this->i == other.i)
&& (this->j == other.j)
&& (this->k == other.k)
&& (this->depth == other.depth)
&& (this->dimensions == other.dimensions)
&& (this->props == other.props)
;
}
Opm::CompletedCells::Cell
Opm::CompletedCells::Cell::serializationTestObject()
{
Cell cell { 0, 1, 1, 1 };
cell.props = Props::serializationTestObject();
cell.depth = 12345;
cell.dimensions = {1.0,2.0,3.0};
return cell;
}
Opm::CompletedCells::CompletedCells(std::size_t nx, std::size_t ny, std::size_t nz)
: dims(nx,ny,nz)
: dims(nx, ny, nz)
{}
Opm::CompletedCells::CompletedCells(const Opm::GridDims& dims_)
:dims(dims_)
: dims(dims_)
{}
const Opm::CompletedCells::Cell& Opm::CompletedCells::get(std::size_t i, std::size_t j, std::size_t k) const {
auto g = this->dims.getGlobalIndex(i,j,k);
return this->cells.at(g);
const Opm::CompletedCells::Cell&
Opm::CompletedCells::get(std::size_t i, std::size_t j, std::size_t k) const
{
return this->cells.at(this->dims.getGlobalIndex(i, j, k));
}
std::pair<bool, Opm::CompletedCells::Cell&>
Opm::CompletedCells::try_get(std::size_t i, std::size_t j, std::size_t k)
{
const auto g = this->dims.getGlobalIndex(i, j, k);
std::pair<bool, Opm::CompletedCells::Cell&> Opm::CompletedCells::try_get(std::size_t i, std::size_t j, std::size_t k) {
auto g = this->dims.getGlobalIndex(i,j,k);
auto iter = this->cells.find(g);
if (iter != this->cells.end())
return {true, iter->second};
const auto& [pos, inserted] = this->cells.try_emplace(g, g, i, j, k);
this->cells.emplace(g, Cell{g,i,j,k});
return {false, this->cells.at(g)};
return { !inserted, pos->second };
}
bool Opm::CompletedCells::operator==(const Opm::CompletedCells& other) const {
return this->dims == other.dims &&
this->cells == other.cells;
bool Opm::CompletedCells::operator==(const Opm::CompletedCells& other) const
{
return (this->dims == other.dims)
&& (this->cells == other.cells)
;
}
Opm::CompletedCells Opm::CompletedCells::serializationTestObject() {
Opm::CompletedCells
Opm::CompletedCells::serializationTestObject()
{
Opm::CompletedCells cells(2,3,4);
cells.cells.emplace(7, Opm::CompletedCells::Cell::serializationTestObject());
return cells;
}
std::size_t Opm::CompletedCells::Cell::active_index() const{
std::size_t Opm::CompletedCells::Cell::active_index() const
{
return this->props.value().active_index;
}
bool Opm::CompletedCells::Cell::is_active() const{
bool Opm::CompletedCells::Cell::is_active() const
{
return this->props.has_value();
}

View File

@ -17,16 +17,24 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <fmt/format.h>
#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
#include <opm/input/eclipse/Schedule/ScheduleGrid.hpp>
#include <opm/input/eclipse/Schedule/CompletedCells.hpp>
#include <opm/input/eclipse/EclipseState/Grid/EclipseGrid.hpp>
#include <opm/input/eclipse/EclipseState/Grid/FieldPropsManager.hpp>
Opm::ScheduleGrid::ScheduleGrid(const Opm::EclipseGrid& ecl_grid, const Opm::FieldPropsManager& fpm, Opm::CompletedCells& completed_cells)
: grid(&ecl_grid)
, fp(&fpm)
, cells(completed_cells)
#include <cstddef>
#include <string>
#include <fmt/format.h>
Opm::ScheduleGrid::ScheduleGrid(const Opm::EclipseGrid& ecl_grid,
const Opm::FieldPropsManager& fpm,
Opm::CompletedCells& completed_cells)
: grid { &ecl_grid }
, fp { &fpm }
, cells { completed_cells }
{}
Opm::ScheduleGrid::ScheduleGrid(Opm::CompletedCells& completed_cells)
@ -34,45 +42,58 @@ Opm::ScheduleGrid::ScheduleGrid(Opm::CompletedCells& completed_cells)
{}
namespace {
double try_get_value(const Opm::FieldPropsManager& fp, const std::string& kw, std::size_t active_index) {
if (fp.has_double(kw))
return fp.try_get<double>(kw)->at(active_index);
else
double try_get_value(const Opm::FieldPropsManager& fp,
const std::string& kw,
const std::size_t active_index)
{
if (! fp.has_double(kw)) {
throw std::logic_error(fmt::format("FieldPropsManager is missing keyword '{}'", kw));
}
double try_get_ntg_value(const Opm::FieldPropsManager& fp, const std::string& kw, std::size_t active_index){
if (fp.has_double(kw))
return fp.try_get<double>(kw)->at(active_index);
else
return 1.0;
}
double try_get_ntg_value(const Opm::FieldPropsManager& fp,
const std::string& kw,
const std::size_t active_index)
{
return fp.has_double(kw)
? fp.try_get<double>(kw)->at(active_index)
: 1.0;
}
}
const Opm::CompletedCells::Cell& Opm::ScheduleGrid::get_cell(std::size_t i, std::size_t j, std::size_t k) const {
if (this->grid) {
auto [valid, cell] = this->cells.try_get(i,j,k);
if (!valid) {
cell.depth = this->grid->getCellDepth(i,j,k);
cell.dimensions = this->grid->getCellDimensions(i,j,k);
if (this->grid->cellActive(i,j,k)){
CompletedCells::Cell::Props props;
const Opm::CompletedCells::Cell&
Opm::ScheduleGrid::get_cell(std::size_t i, std::size_t j, std::size_t k) const
{
if (this->grid == nullptr) {
return this->cells.get(i, j, k);
}
props.active_index = this->grid->getActiveIndex(i,j,k);
auto [valid, cellRef] = this->cells.try_get(i, j, k);
if (!valid) {
cellRef.depth = this->grid->getCellDepth(i, j, k);
cellRef.dimensions = this->grid->getCellDimensions(i, j, k);
if (this->grid->cellActive(i, j, k)) {
auto& props = cellRef.props.emplace(CompletedCells::Cell::Props{});
props.active_index = this->grid->getActiveIndex(i, j, k);
props.permx = try_get_value(*this->fp, "PERMX", props.active_index);
props.permy = try_get_value(*this->fp, "PERMY", props.active_index);
props.permz = try_get_value(*this->fp, "PERMZ", props.active_index);
props.poro = try_get_value(*this->fp, "PORO", props.active_index);
props.satnum = this->fp->get_int("SATNUM").at(props.active_index);
props.pvtnum = this->fp->get_int("PVTNUM").at(props.active_index);
props.ntg = try_get_ntg_value(*this->fp, "NTG", props.active_index);
cell.props = props;
}
}
return cell;
} else
return this->cells.get(i,j,k);
return cellRef;
}
const Opm::EclipseGrid* Opm::ScheduleGrid::get_grid() const {
const Opm::EclipseGrid* Opm::ScheduleGrid::get_grid() const
{
return this->grid;
}

View File

@ -263,6 +263,10 @@ COPY
'PERMX' 'PERMZ' /
'PERMX' 'PERMY' /
/
PORO
1000*0.3 /
SCHEDULE
COMPDAT
@ -288,6 +292,9 @@ COPY
'PERMX' 'PERMY' /
/
PORO
1000*0.3 /
SCHEDULE
COMPDAT

View File

@ -14,6 +14,8 @@ COPY
PERMX PERMZ /
/
PORO
1000*0.3 /
DX
1000*1 /

View File

@ -16,6 +16,9 @@ COPY
PERMX PERMZ /
/
PORO
48000*0.3 /
SCHEDULE
WELSPECS

View File

@ -17,6 +17,9 @@ COPY
PERMX PERMZ /
/
PORO
27000*0.3 /
SCHEDULE

View File

@ -15,6 +15,9 @@ COPY
PERMX PERMZ /
/
PORO
27000*0.3 /
SCHEDULE

View File

@ -17,6 +17,8 @@ COPY
PERMX PERMZ /
/
PORO
800*0.3 /
SCHEDULE

View File

@ -16,6 +16,8 @@ COPY
PERMX PERMZ /
/
PORO
27000*0.3 /
SCHEDULE

View File

@ -13,6 +13,9 @@ COPY
PERMX PERMZ /
/
PORO
72000*0.3 /
SCHEDULE
-- 0