Track AQUNUM Records for Restart Purposes

When forming restart arrays for numerical aquifers (IAQN and RAQN)
we need to know the total number of records in the input AQUNUM
keyword as well as the record ID of each individual record.  This
commit adds a tracking mechanism for this information.

While here, also add unit tests that verify that we've correctly
interpreted the AQUNUM records.
This commit is contained in:
Bård Skaflestad 2021-05-18 12:30:19 +02:00
parent f58a314f6a
commit 8eb9bee6ed
5 changed files with 116 additions and 37 deletions

View File

@ -20,6 +20,7 @@
#ifndef OPM_NUMERICALAQUIFERCELL_HPP
#define OPM_NUMERICALAQUIFERCELL_HPP
#include <cstddef>
#include <optional>
namespace Opm {
@ -28,19 +29,20 @@ namespace Opm {
class FieldPropsManager;
struct NumericalAquiferCell {
NumericalAquiferCell(const DeckRecord&, const EclipseGrid&, const FieldPropsManager&);
NumericalAquiferCell(const std::size_t record_id_, const DeckRecord&, const EclipseGrid&, const FieldPropsManager&);
NumericalAquiferCell() = default;
size_t aquifer_id;
size_t I, J, K;
double area;
double length;
double porosity;
double permeability;
double depth;
std::optional<double> init_pressure;
int pvttable;
int sattable;
size_t global_index;
std::size_t aquifer_id{};
std::size_t I{}, J{}, K{};
double area{};
double length{};
double porosity{};
double permeability{};
double depth{};
std::optional<double> init_pressure{};
int pvttable{};
int sattable{};
std::size_t global_index{};
std::size_t record_id{};
double cellVolume() const;
double poreVolume() const;
@ -62,6 +64,7 @@ namespace Opm {
serializer(this->pvttable);
serializer(this->sattable);
serializer(this->global_index);
serializer(this->record_id);
}
};
}

View File

@ -20,8 +20,10 @@
#ifndef OPM_NUMERICALAQUIFERS_HPP
#define OPM_NUMERICALAQUIFERS_HPP
#include <unordered_map>
#include <map>
#include <unordered_map>
#include <stddef.h>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/SingleNumericalAquifer.hpp>
@ -36,6 +38,7 @@ namespace Opm {
NumericalAquifers() = default;
NumericalAquifers(const Deck& deck, const EclipseGrid& grid, const FieldPropsManager& field_props);
int numRecords() const { return static_cast<int>(this->m_num_records); }
size_t size() const;
bool hasAquifer(size_t aquifer_id) const;
const SingleNumericalAquifer& getAquifer(size_t aquifer_id) const;
@ -59,10 +62,12 @@ namespace Opm {
void serializeOp(Serializer& serializer)
{
serializer.map(this->m_aquifers);
serializer(this->m_num_records);
}
private:
std::map<size_t, SingleNumericalAquifer> m_aquifers;
std::map<size_t, SingleNumericalAquifer> m_aquifers{};
size_t m_num_records{0};
void addAquiferCell(const NumericalAquiferCell& aqu_cell);
};

View File

@ -27,21 +27,24 @@
namespace Opm {
using AQUNUM = ParserKeywords::AQUNUM;
NumericalAquiferCell::NumericalAquiferCell(const DeckRecord& record, const EclipseGrid& grid, const FieldPropsManager& field_props)
: aquifer_id( record.getItem<AQUNUM::AQUIFER_ID>().get<int>(0) )
, I ( record.getItem<AQUNUM::I>().get<int>(0) - 1 )
, J ( record.getItem<AQUNUM::J>().get<int>(0) - 1 )
, K ( record.getItem<AQUNUM::K>().get<int>(0) - 1 )
, area (record.getItem<AQUNUM::CROSS_SECTION>().getSIDouble(0) )
, length ( record.getItem<AQUNUM::LENGTH>().getSIDouble(0) )
, permeability( record.getItem<AQUNUM::PERM>().getSIDouble(0) )
NumericalAquiferCell::NumericalAquiferCell(const std::size_t record_id_,
const DeckRecord& record,
const EclipseGrid& grid,
const FieldPropsManager& field_props)
: aquifer_id( record.getItem<AQUNUM::AQUIFER_ID>().get<int>(0) )
, I ( record.getItem<AQUNUM::I>().get<int>(0) - 1 )
, J ( record.getItem<AQUNUM::J>().get<int>(0) - 1 )
, K ( record.getItem<AQUNUM::K>().get<int>(0) - 1 )
, area (record.getItem<AQUNUM::CROSS_SECTION>().getSIDouble(0) )
, length ( record.getItem<AQUNUM::LENGTH>().getSIDouble(0) )
, permeability( record.getItem<AQUNUM::PERM>().getSIDouble(0) )
{
const auto& poro = field_props.get_double("PORO");
const auto& pvtnum = field_props.get_int("PVTNUM");
const auto& satnum = field_props.get_int("SATNUM");
this->global_index = grid.getGlobalIndex(I, J, K);
const size_t active_index = grid.activeIndex(this->global_index);
const std::size_t active_index = grid.activeIndex(this->global_index);
if ( !record.getItem<AQUNUM::PORO>().defaultApplied(0) ) {
this->porosity = record.getItem<AQUNUM::PORO>().getSIDouble(0);
@ -70,6 +73,8 @@ namespace Opm {
} else {
this->sattable = satnum[active_index];
}
this->record_id = record_id_;
}
double NumericalAquiferCell::cellVolume() const {
@ -89,7 +94,8 @@ namespace Opm {
this->init_pressure == other.init_pressure &&
this->pvttable == other.pvttable &&
this->sattable == other.sattable &&
this->global_index == other.global_index;
this->global_index == other.global_index &&
this->record_id == other.record_id;
}
double NumericalAquiferCell::poreVolume() const {

View File

@ -17,38 +17,42 @@
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include <opm/parser/eclipse/Parser/ParserKeywords/A.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquifers.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquiferCell.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/SingleNumericalAquifer.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
#include <fmt/format.h>
#include <opm/common/utility/OpmInputError.hpp>
#include <opm/common/OpmLog/OpmLog.hpp>
#include <opm/common/OpmLog/KeywordLocation.hpp>
#include <opm/parser/eclipse/Deck/Deck.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
#include <opm/parser/eclipse/EclipseState/Grid/EclipseGrid.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquiferCell.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/NumericalAquifers.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/NumericalAquifer/SingleNumericalAquifer.hpp>
#include <opm/parser/eclipse/Parser/ParserKeywords/A.hpp>
#include <set>
#include <fmt/format.h>
#include <cstddef>
#include <unordered_set>
namespace Opm {
NumericalAquifers::NumericalAquifers(const Deck& deck, const EclipseGrid& grid,
const FieldPropsManager& field_props) {
const FieldPropsManager& field_props)
{
using AQUNUM=ParserKeywords::AQUNUM;
if ( !deck.hasKeyword<AQUNUM>() ) return;
std::set<size_t> cells;
std::unordered_set<std::size_t> cells;
// there might be multiple keywords of keyword AQUNUM, it is not totally
// clear about the rules here. For now, we take care of all the keywords
const auto& aqunum_keywords = deck.getKeywordList<AQUNUM>();
for (const auto& keyword : aqunum_keywords) {
for (const auto& record : *keyword) {
const NumericalAquiferCell aqu_cell(record, grid, field_props);
const NumericalAquiferCell aqu_cell(this->m_num_records++, record, grid, field_props);
if (cells.count(aqu_cell.global_index) > 0) {
auto error = fmt::format("Numerical aquifer cell at ({}, {}, {}) is declared more than once",
aqu_cell.I + 1, aqu_cell.J + 1, aqu_cell.K + 1);
@ -115,7 +119,8 @@ namespace Opm {
}
bool NumericalAquifers::operator==(const NumericalAquifers& other) const {
return this->m_aquifers == other.m_aquifers;
return (this->m_aquifers == other.m_aquifers)
&& (this->m_num_records == other.m_num_records);
}
size_t NumericalAquifers::size() const {
@ -196,4 +201,4 @@ namespace Opm {
aquifer.postProcessConnections(grid, actnum);
}
}
}
}

View File

@ -25,6 +25,8 @@ along with OPM. If not, see <http://www.gnu.org/licenses/>.
#include <opm/parser/eclipse/EclipseState/Aquifer/Aquifetp.hpp>
#include <opm/parser/eclipse/EclipseState/Aquifer/AquiferConfig.hpp>
#include <opm/parser/eclipse/Units/Units.hpp>
using namespace Opm;
@ -717,6 +719,64 @@ BOOST_AUTO_TEST_CASE(NumericalAquiferTest){
const Opm::EclipseGrid& grid = ecl_state.getInputGrid();
Opm::NumericalAquifers num_aqu{numaquifer_deck, grid, ecl_state.fieldProps()};
BOOST_CHECK_EQUAL(num_aqu.numRecords(), 3);
{
const auto mD = unit::convert::from(1.0, prefix::milli*unit::darcy);
const auto* c1 = num_aqu.getAquifer(1).getCellPrt(0);
const auto* c2 = num_aqu.getAquifer(1).getCellPrt(1);
const auto* c3 = num_aqu.getAquifer(1).getCellPrt(2);
BOOST_CHECK_EQUAL(c1->record_id, std::size_t{0});
BOOST_CHECK_EQUAL(c1->I, std::size_t{0});
BOOST_CHECK_EQUAL(c1->J, std::size_t{0});
BOOST_CHECK_EQUAL(c1->K, std::size_t{0});
BOOST_CHECK_CLOSE(c1->area, 1.0e6, 1.0e-10);
BOOST_CHECK_CLOSE(c1->length, 10.0e3, 1.0e-10);
BOOST_CHECK_CLOSE(c1->porosity, 0.25, 1.0e-10);
BOOST_CHECK_CLOSE(c1->permeability, 400*mD, 1.0e-10);
BOOST_CHECK_CLOSE(c1->depth, 2585.0, 1.0e-10);
BOOST_CHECK_MESSAGE(c1->init_pressure.has_value(), "Cell 1 must have an initial pressure");
BOOST_CHECK_CLOSE(c1->init_pressure.value(), 285.0*unit::barsa, 1.0e-10);
BOOST_CHECK_EQUAL(c1->pvttable, 2);
BOOST_CHECK_EQUAL(c1->sattable, 2);
BOOST_CHECK_EQUAL(c1->global_index, 0);
BOOST_CHECK_EQUAL(c2->record_id, std::size_t{1});
BOOST_CHECK_EQUAL(c2->I, std::size_t{2});
BOOST_CHECK_EQUAL(c2->J, std::size_t{0});
BOOST_CHECK_EQUAL(c2->K, std::size_t{0});
BOOST_CHECK_CLOSE(c2->area, 1.5e6, 1.0e-10);
BOOST_CHECK_CLOSE(c2->length, 20.0e3, 1.0e-10);
BOOST_CHECK_CLOSE(c2->porosity, 0.24, 1.0e-10);
BOOST_CHECK_CLOSE(c2->permeability, 600*mD, 1.0e-10);
BOOST_CHECK_CLOSE(c2->depth, 2585.0, 1.0e-10);
BOOST_CHECK_MESSAGE(c2->init_pressure.has_value(), "Cell 2 must have an initial pressure");
BOOST_CHECK_CLOSE(c2->init_pressure.value(), 285.0*unit::barsa, 1.0e-10);
BOOST_CHECK_EQUAL(c2->pvttable, 3);
BOOST_CHECK_EQUAL(c2->sattable, 1);
BOOST_CHECK_EQUAL(c2->global_index, 2);
BOOST_CHECK_EQUAL(c3->record_id, std::size_t{2});
BOOST_CHECK_EQUAL(c3->I, std::size_t{3});
BOOST_CHECK_EQUAL(c3->J, std::size_t{0});
BOOST_CHECK_EQUAL(c3->K, std::size_t{0});
BOOST_CHECK_CLOSE(c3->area, 2.0e6, 1.0e-10);
BOOST_CHECK_CLOSE(c3->length, 30.0e3, 1.0e-10);
BOOST_CHECK_CLOSE(c3->porosity, 0.25, 1.0e-10);
BOOST_CHECK_CLOSE(c3->permeability, 700*mD, 1.0e-10);
BOOST_CHECK_CLOSE(c3->depth, 2585.0, 1.0e-10);
BOOST_CHECK_MESSAGE(c3->init_pressure.has_value(), "Cell 3 must have an initial pressure");
BOOST_CHECK_CLOSE(c3->init_pressure.value(), 285.0*unit::barsa, 1.0e-10);
BOOST_CHECK_EQUAL(c3->pvttable, 1);
BOOST_CHECK_EQUAL(c3->sattable, 3);
BOOST_CHECK_EQUAL(c3->global_index, 3);
}
// using processed actnum for numerical aquifer connection generation
std::vector<int> new_actnum(360, 1);
new_actnum[0] = 0;