changed: avoid typetag dependence in equil initializer

This commit is contained in:
Arne Morten Kvarving 2022-08-09 12:16:27 +02:00
parent 624d7e51cd
commit b399f72777
3 changed files with 79 additions and 51 deletions

View File

@ -100,12 +100,21 @@ public:
const auto& eclState = vanguard.eclState(); const auto& eclState = vanguard.eclState();
unsigned numElems = vanguard.grid().size(0); unsigned numElems = vanguard.grid().size(0);
using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
using Grid = GetPropType<TypeTag, Properties::Grid>;
using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
using Initializer = EQUIL::DeckDependent::InitialStateComputer<FluidSystem,
Grid,
GridView,
ElementMapper,
CartesianIndexMapper>;
EQUIL::DeckDependent::InitialStateComputer<TypeTag> initialState(materialLawManager, Initializer initialState(materialLawManager,
eclState, eclState,
vanguard, vanguard.grid(),
vanguard.cartesianMapper(), vanguard.gridView(),
simulator.problem().gravity()[dimWorld - 1]); vanguard.cartesianMapper(),
simulator.problem().gravity()[dimWorld - 1]);
// copy the result into the array of initial fluid states // copy the result into the array of initial fluid states
initialFluidStates_.resize(numElems); initialFluidStates_.resize(numElems);

View File

@ -1581,34 +1581,32 @@ equilnum(const EclipseState& eclipseState,
return eqlnum; return eqlnum;
} }
template<class TypeTag> template<class FluidSystem,
class Grid,
class GridView,
class ElementMapper,
class CartesianIndexMapper>
class InitialStateComputer class InitialStateComputer
{ {
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
using GridView = GetPropType<TypeTag, Properties::GridView>;
using Element = typename GridView::template Codim<0>::Entity; using Element = typename GridView::template Codim<0>::Entity;
using ElementMapper = GetPropType<TypeTag, Properties::ElementMapper>;
using Vanguard = GetPropType<TypeTag, Properties::Vanguard>;
using Grid = GetPropType<TypeTag, Properties::Grid>;
using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
public: public:
template<class MaterialLawManager> template<class MaterialLawManager>
InitialStateComputer(MaterialLawManager& materialLawManager, InitialStateComputer(MaterialLawManager& materialLawManager,
const Opm::EclipseState& eclipseState, const Opm::EclipseState& eclipseState,
const Vanguard& vanguard, const Grid& grid,
const GridView& gridView,
const CartesianIndexMapper& cartMapper, const CartesianIndexMapper& cartMapper,
const double grav = Opm::unit::gravity, const double grav = Opm::unit::gravity,
const bool applySwatInit = true) const bool applySwatInit = true)
: temperature_(vanguard.grid().size(/*codim=*/0)), : temperature_(grid.size(/*codim=*/0)),
saltConcentration_(vanguard.grid().size(/*codim=*/0)), saltConcentration_(grid.size(/*codim=*/0)),
saltSaturation_(vanguard.grid().size(/*codim=*/0)), saltSaturation_(grid.size(/*codim=*/0)),
pp_(FluidSystem::numPhases, pp_(FluidSystem::numPhases,
std::vector<double>(vanguard.grid().size(/*codim=*/0))), std::vector<double>(grid.size(/*codim=*/0))),
sat_(FluidSystem::numPhases, sat_(FluidSystem::numPhases,
std::vector<double>(vanguard.grid().size(/*codim=*/0))), std::vector<double>(grid.size(/*codim=*/0))),
rs_(vanguard.grid().size(/*codim=*/0)), rs_(grid.size(/*codim=*/0)),
rv_(vanguard.grid().size(/*codim=*/0)), rv_(grid.size(/*codim=*/0)),
cartesianIndexMapper_(cartMapper) cartesianIndexMapper_(cartMapper)
{ {
//Check for presence of kw SWATINIT //Check for presence of kw SWATINIT
@ -1621,8 +1619,7 @@ public:
// Querry cell depth, cell top-bottom. // Querry cell depth, cell top-bottom.
// numerical aquifer cells might be specified with different depths. // numerical aquifer cells might be specified with different depths.
const auto& num_aquifers = eclipseState.aquifer().numericalAquifers(); const auto& num_aquifers = eclipseState.aquifer().numericalAquifers();
updateCellProps_(vanguard, num_aquifers); updateCellProps_(gridView, num_aquifers);
const Grid& grid = vanguard.grid();
// Get the equilibration records. // Get the equilibration records.
const std::vector<EquilRecord> rec = getEquil(eclipseState); const std::vector<EquilRecord> rec = getEquil(eclipseState);
@ -1743,7 +1740,7 @@ public:
calcPressSatRsRv(eqlmap, rec, materialLawManager, comm, grav); calcPressSatRsRv(eqlmap, rec, materialLawManager, comm, grav);
// modify the pressure and saturation for numerical aquifer cells // modify the pressure and saturation for numerical aquifer cells
applyNumericalAquifers_(vanguard, num_aquifers, eclipseState.runspec().co2Storage()); applyNumericalAquifers_(gridView, num_aquifers, eclipseState.runspec().co2Storage());
// Modify oil pressure in no-oil regions so that the pressures of present phases can // Modify oil pressure in no-oil regions so that the pressures of present phases can
// be recovered from the oil pressure and capillary relations. // be recovered from the oil pressure and capillary relations.
@ -1837,10 +1834,9 @@ private:
std::vector<std::pair<double,double>> cellZSpan_; std::vector<std::pair<double,double>> cellZSpan_;
std::vector<std::pair<double,double>> cellZMinMax_; std::vector<std::pair<double,double>> cellZMinMax_;
void updateCellProps_(const Vanguard& vanguard, void updateCellProps_(const GridView& gridView,
const NumericalAquifers& aquifer) const NumericalAquifers& aquifer)
{ {
const auto& gridView = vanguard.gridView();
ElementMapper elemMapper(gridView, Dune::mcmgElementLayout()); ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
//const auto& elemMapper = gridView.elementMapper(); //const auto& elemMapper = gridView.elementMapper();
int numElements = gridView.size(/*codim=*/0); int numElements = gridView.size(/*codim=*/0);
@ -1873,13 +1869,12 @@ private:
} }
} }
void applyNumericalAquifers_(const Vanguard& vanguard, void applyNumericalAquifers_(const GridView& gridView,
const NumericalAquifers& aquifer, const NumericalAquifers& aquifer,
const bool co2store) const bool co2store)
{ {
const auto num_aqu_cells = aquifer.allAquiferCells(); const auto num_aqu_cells = aquifer.allAquiferCells();
if (num_aqu_cells.empty()) return; if (num_aqu_cells.empty()) return;
const auto& gridView = vanguard.gridView();
ElementMapper elemMapper(gridView, Dune::mcmgElementLayout()); ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
//const auto& elemMapper = vanguard.gridView().elementMapper();; //const auto& elemMapper = vanguard.gridView().elementMapper();;
auto elemIt = gridView.template begin</*codim=*/0>(); auto elemIt = gridView.template begin</*codim=*/0>();

View File

@ -220,7 +220,6 @@ struct EquilFixture {
Dune::MPIHelper::instance(argc, argv); Dune::MPIHelper::instance(argc, argv);
#endif #endif
Opm::EclGenericVanguard::setCommunication(std::make_unique<Opm::Parallel::Communication>()); Opm::EclGenericVanguard::setCommunication(std::make_unique<Opm::Parallel::Communication>());
using TypeTag = Opm::Properties::TTag::TestEquilTypeTag;
Opm::BlackoilModelParametersEbos<TypeTag>::registerParameters(); Opm::BlackoilModelParametersEbos<TypeTag>::registerParameters();
Opm::Parameters::registerParam<TypeTag, bool>("EnableTerminalOutput", Opm::Parameters::registerParam<TypeTag, bool>("EnableTerminalOutput",
"EnableTerminalOutput", "EnableTerminalOutput",
@ -228,6 +227,18 @@ struct EquilFixture {
"Dummy added for the well model to compile."); "Dummy added for the well model to compile.");
Opm::registerAllParameters_<TypeTag>(); Opm::registerAllParameters_<TypeTag>();
} }
using TypeTag = Opm::Properties::TTag::TestEquilTypeTag;
using FluidSystem = Opm::GetPropType<TypeTag, Opm::Properties::FluidSystem>;
using Grid = Opm::GetPropType<TypeTag, Opm::Properties::Grid>;
using GridView = Opm::GetPropType<TypeTag, Opm::Properties::GridView>;
using ElementMapper = Opm::GetPropType<TypeTag, Opm::Properties::ElementMapper>;
using CartesianIndexMapper = Dune::CartesianIndexMapper<Grid>;
using Initializer = Opm::EQUIL::DeckDependent::InitialStateComputer<FluidSystem,
Grid,
GridView,
ElementMapper,
CartesianIndexMapper>;
}; };
} }
@ -501,10 +512,11 @@ BOOST_AUTO_TEST_CASE(DeckAllDead)
const auto& eclipseState = simulator->vanguard().eclState(); const auto& eclipseState = simulator->vanguard().eclState();
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), eclipseState,
eclipseState, simulator->vanguard(), simulator->vanguard().grid(),
simulator->vanguard().cartesianMapper(), 10.0); simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 10.0);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);
@ -581,9 +593,11 @@ BOOST_AUTO_TEST_CASE(DeckWithCapillary)
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 10.0); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 10.0);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
@ -621,9 +635,11 @@ BOOST_AUTO_TEST_CASE(DeckWithCapillaryOverlap)
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 9.80665); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 9.80665);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);
@ -682,9 +698,11 @@ BOOST_AUTO_TEST_CASE(DeckWithLiveOil)
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
// Initialize the fluid system // Initialize the fluid system
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 9.80665); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 9.80665);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);
@ -759,9 +777,11 @@ BOOST_AUTO_TEST_CASE(DeckWithLiveGas)
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 9.80665); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 9.80665);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);
@ -839,9 +859,11 @@ BOOST_AUTO_TEST_CASE(DeckWithRSVDAndRVVD)
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 9.80665); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 9.80665);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);
@ -939,9 +961,11 @@ BOOST_AUTO_TEST_CASE(DeckWithPBVDAndPDVD)
Opm::GridManager gm(eclipseState.getInputGrid()); Opm::GridManager gm(eclipseState.getInputGrid());
const UnstructuredGrid& grid = *(gm.c_grid()); const UnstructuredGrid& grid = *(gm.c_grid());
Opm::EQUIL::DeckDependent::InitialStateComputer<TypeTag> comp(*simulator->problem().materialLawManager(), EquilFixture::Initializer comp(*simulator->problem().materialLawManager(),
eclipseState, simulator->vanguard(), eclipseState,
simulator->vanguard().cartesianMapper(), 9.80665); simulator->vanguard().grid(),
simulator->vanguard().gridView(),
simulator->vanguard().cartesianMapper(), 9.80665);
const auto& pressures = comp.press(); const auto& pressures = comp.press();
BOOST_REQUIRE_EQUAL(pressures.size(), 3U); BOOST_REQUIRE_EQUAL(pressures.size(), 3U);
BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells); BOOST_REQUIRE_EQUAL(int(pressures[0].size()), grid.number_of_cells);