mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-11 23:05:33 -06:00
Add trivial saltvd table for cases without saltvd given in the deck
This commit is contained in:
parent
bd9186b41d
commit
491b532a6a
@ -178,7 +178,8 @@ private:
|
||||
density(const double depth,
|
||||
const double press) const
|
||||
{
|
||||
double saltConcentration = saltVdTable_.eval(depth);
|
||||
// The initializing algorithm can give depths outside the range due to numerical noise i.e. we extrapolate
|
||||
double saltConcentration = saltVdTable_.eval(depth, /*extrapolate=*/true);
|
||||
double rho = FluidSystem::waterPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp_, press, saltConcentration);
|
||||
rho *= FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
|
||||
return rho;
|
||||
@ -1736,17 +1737,28 @@ private:
|
||||
template <class RMap>
|
||||
void updateInitialSaltConcentration_(const Opm::EclipseState& eclState, const RMap& reg, const Grid& grid)
|
||||
{
|
||||
const int numEquilReg = rsFunc_.size();
|
||||
saltVdTable_.resize(numEquilReg);
|
||||
const auto& tables = eclState.getTableManager();
|
||||
const Opm::TableContainer& saltvdTables = tables.getSaltvdTables();
|
||||
saltVdTable_.resize(saltvdTables.size());
|
||||
for (size_t i = 0; i < saltvdTables.size(); ++i) {
|
||||
const Opm::SaltvdTable& saltvdTable = saltvdTables.getTable<Opm::SaltvdTable>(i);
|
||||
saltVdTable_[i].setXYContainers(saltvdTable.getDepthColumn(), saltvdTable.getSaltColumn());
|
||||
|
||||
const auto& cells = reg.cells(i);
|
||||
for (const auto& cell : cells) {
|
||||
const double depth = UgGridHelpers::cellCenterDepth(grid, cell);
|
||||
this->saltConcentration_[cell] = saltVdTable_[i].eval(depth);
|
||||
// If no saltvd table is given, we create a trivial table for the density calculations
|
||||
if (saltvdTables.empty()) {
|
||||
std::vector<double> x = {0.0,1.0};
|
||||
std::vector<double> y = {0.0,0.0};
|
||||
for (size_t i = 0; i < numEquilReg; ++i) {
|
||||
saltVdTable_[i].setXYContainers(x,y);
|
||||
}
|
||||
} else {
|
||||
for (size_t i = 0; i < saltvdTables.size(); ++i) {
|
||||
const Opm::SaltvdTable& saltvdTable = saltvdTables.getTable<Opm::SaltvdTable>(i);
|
||||
saltVdTable_[i].setXYContainers(saltvdTable.getDepthColumn(), saltvdTable.getSaltColumn());
|
||||
|
||||
const auto& cells = reg.cells(i);
|
||||
for (const auto& cell : cells) {
|
||||
const double depth = UgGridHelpers::cellCenterDepth(grid, cell);
|
||||
this->saltConcentration_[cell] = saltVdTable_[i].eval(depth);
|
||||
}
|
||||
}
|
||||
}
|
||||
}
|
||||
|
@ -182,6 +182,11 @@ void test_PhasePressure()
|
||||
using TypeTag = TTAG(TestEquilTypeTag);
|
||||
using FluidSystem = GET_PROP_TYPE(TypeTag, FluidSystem);
|
||||
|
||||
using TabulatedFunction = typename Opm::Tabulated1DFunction<double>;
|
||||
std::vector<double> x = {0.0,100.0};
|
||||
std::vector<double> y = {0.0,0.0};
|
||||
TabulatedFunction trivialSaltVdTable = TabulatedFunction(2, x, y);
|
||||
|
||||
auto simulator = initSimulator<TypeTag>("equil_base.DATA");
|
||||
initDefaultFluidSystem<TypeTag>();
|
||||
|
||||
@ -189,6 +194,7 @@ void test_PhasePressure()
|
||||
record,
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0
|
||||
};
|
||||
|
||||
@ -235,26 +241,35 @@ void test_CellSubset()
|
||||
const Opm::EquilRecord record[] = { mkEquilRecord( 0, 1e5, 2.5, -0.075e5, 0, 0 ),
|
||||
mkEquilRecord( 5, 1.35e5, 7.5, -0.225e5, 5, 0 ) };
|
||||
|
||||
using TabulatedFunction = typename Opm::Tabulated1DFunction<double>;
|
||||
std::vector<double> x = {0.0,100.0};
|
||||
std::vector<double> y = {0.0,0.0};
|
||||
TabulatedFunction trivialSaltVdTable = TabulatedFunction(2, x, y);
|
||||
|
||||
const Opm::EQUIL::EquilReg region[] =
|
||||
{
|
||||
Opm::EQUIL::EquilReg(record[0],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[0],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[1],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[1],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
};
|
||||
|
||||
@ -330,26 +345,35 @@ void test_RegMapping()
|
||||
auto simulator = initSimulator<TypeTag>("equil_base.DATA");
|
||||
initDefaultFluidSystem<TypeTag>();
|
||||
|
||||
using TabulatedFunction = typename Opm::Tabulated1DFunction<double>;
|
||||
std::vector<double> x = {0.0,100.0};
|
||||
std::vector<double> y = {0.0,0.0};
|
||||
TabulatedFunction trivialSaltVdTable = TabulatedFunction(2, x, y);
|
||||
|
||||
const Opm::EQUIL::EquilReg region[] =
|
||||
{
|
||||
Opm::EQUIL::EquilReg(record[0],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[0],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[1],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
,
|
||||
Opm::EQUIL::EquilReg(record[1],
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
std::make_shared<Opm::EQUIL::Miscibility::NoMixing>(),
|
||||
trivialSaltVdTable,
|
||||
0)
|
||||
};
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user