2018-10-16 08:59:16 -05:00
|
|
|
#include <opm/grid/utility/cartesianToCompressed.hpp>
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
BlackoilAquiferModel<TypeTag>::
|
|
|
|
BlackoilAquiferModel(Simulator& simulator)
|
|
|
|
: simulator_(simulator)
|
|
|
|
{
|
|
|
|
init();
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>::initialSolutionApplied()
|
|
|
|
{
|
|
|
|
if(aquiferCarterTracyActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_CarterTracy.begin(); aquifer != aquifers_CarterTracy.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->initialSolutionApplied();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(aquiferFetkovichActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_Fetkovich.begin(); aquifer != aquifers_Fetkovich.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->initialSolutionApplied();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>::beginEpisode()
|
|
|
|
{ }
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>::beginIteration()
|
|
|
|
{ }
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void BlackoilAquiferModel<TypeTag>:: beginTimeStep()
|
|
|
|
{
|
|
|
|
if(aquiferCarterTracyActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_CarterTracy.begin(); aquifer != aquifers_CarterTracy.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->beginTimeStep();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(aquiferFetkovichActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_Fetkovich.begin(); aquifer != aquifers_Fetkovich.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->beginTimeStep();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
template<class Context>
|
|
|
|
void BlackoilAquiferModel<TypeTag>:: addToSource(RateVector& rates, const Context& context, unsigned spaceIdx, unsigned timeIdx) const
|
|
|
|
{
|
|
|
|
if(aquiferCarterTracyActive())
|
|
|
|
{
|
|
|
|
for (auto& aquifer : aquifers_CarterTracy)
|
|
|
|
{
|
|
|
|
aquifer.addToSource(rates, context, spaceIdx, timeIdx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(aquiferFetkovichActive())
|
|
|
|
{
|
|
|
|
for (auto& aquifer : aquifers_Fetkovich)
|
|
|
|
{
|
|
|
|
aquifer.addToSource(rates, context, spaceIdx, timeIdx);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>::endIteration()
|
|
|
|
{ }
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
void BlackoilAquiferModel<TypeTag>:: endTimeStep()
|
|
|
|
{
|
|
|
|
if(aquiferCarterTracyActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_CarterTracy.begin(); aquifer != aquifers_CarterTracy.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->endTimeStep();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(aquiferFetkovichActive())
|
|
|
|
{
|
|
|
|
for (auto aquifer = aquifers_Fetkovich.begin(); aquifer != aquifers_Fetkovich.end(); ++aquifer)
|
|
|
|
{
|
|
|
|
aquifer->endTimeStep();
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>::endEpisode()
|
|
|
|
{ }
|
|
|
|
|
2019-02-25 05:14:38 -06:00
|
|
|
template <typename TypeTag>
|
|
|
|
template <class Restarter>
|
|
|
|
void
|
2019-04-24 10:18:44 -05:00
|
|
|
BlackoilAquiferModel<TypeTag>::serialize(Restarter& /* res */)
|
2019-02-25 05:14:38 -06:00
|
|
|
{
|
|
|
|
// TODO (?)
|
|
|
|
throw std::logic_error("BlackoilAquiferModel::serialize() is not yet implemented");
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename TypeTag>
|
|
|
|
template <class Restarter>
|
|
|
|
void
|
2019-04-24 10:18:44 -05:00
|
|
|
BlackoilAquiferModel<TypeTag>::deserialize(Restarter& /* res */)
|
2019-02-25 05:14:38 -06:00
|
|
|
{
|
|
|
|
// TODO (?)
|
|
|
|
throw std::logic_error("BlackoilAquiferModel::deserialize() is not yet implemented");
|
|
|
|
}
|
|
|
|
|
2018-10-16 08:59:16 -05:00
|
|
|
// Initialize the aquifers in the deck
|
|
|
|
template<typename TypeTag>
|
|
|
|
void
|
|
|
|
BlackoilAquiferModel<TypeTag>:: init()
|
|
|
|
{
|
|
|
|
const auto& deck = this->simulator_.vanguard().deck();
|
|
|
|
if (deck.hasKeyword("AQUCT")) {
|
|
|
|
//updateConnectionIntensiveQuantities();
|
|
|
|
const auto& eclState = this->simulator_.vanguard().eclState();
|
|
|
|
|
|
|
|
// Get all the carter tracy aquifer properties data and put it in aquifers vector
|
|
|
|
const AquiferCT aquiferct = AquiferCT(eclState,deck);
|
|
|
|
const Aquancon aquifer_connect = Aquancon(eclState.getInputGrid(), deck);
|
|
|
|
|
|
|
|
std::vector<AquiferCT::AQUCT_data> aquifersData = aquiferct.getAquifers();
|
|
|
|
std::vector<Aquancon::AquanconOutput> aquifer_connection = aquifer_connect.getAquOutput();
|
|
|
|
|
|
|
|
assert( aquifersData.size() == aquifer_connection.size() );
|
2019-01-03 06:18:01 -06:00
|
|
|
const auto& ugrid = simulator_.vanguard().grid();
|
|
|
|
const auto& gridView = simulator_.gridView();
|
|
|
|
const int number_of_cells = gridView.size(0);
|
|
|
|
|
2018-10-16 08:59:16 -05:00
|
|
|
cartesian_to_compressed_ = cartesianToCompressed(number_of_cells,
|
|
|
|
Opm::UgGridHelpers::globalCell(ugrid));
|
|
|
|
|
|
|
|
for (size_t i = 0; i < aquifersData.size(); ++i)
|
|
|
|
{
|
|
|
|
aquifers_CarterTracy.push_back(
|
2019-01-31 08:26:09 -06:00
|
|
|
AquiferCarterTracy<TypeTag> (aquifer_connection.at(i), cartesian_to_compressed_, this->simulator_ , aquifersData.at(i))
|
2018-10-16 08:59:16 -05:00
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
if(deck.hasKeyword("AQUFETP"))
|
|
|
|
{
|
|
|
|
//updateConnectionIntensiveQuantities();
|
|
|
|
const auto& eclState = this->simulator_.vanguard().eclState();
|
|
|
|
|
|
|
|
// Get all the carter tracy aquifer properties data and put it in aquifers vector
|
|
|
|
const Aquifetp aquifetp = Aquifetp(deck);
|
|
|
|
const Aquancon aquifer_connect = Aquancon(eclState.getInputGrid(), deck);
|
|
|
|
|
|
|
|
std::vector<Aquifetp::AQUFETP_data> aquifersData = aquifetp.getAquifers();
|
|
|
|
std::vector<Aquancon::AquanconOutput> aquifer_connection = aquifer_connect.getAquOutput();
|
|
|
|
|
|
|
|
assert( aquifersData.size() == aquifer_connection.size() );
|
2019-01-03 06:50:18 -06:00
|
|
|
const auto& ugrid = simulator_.vanguard().grid();
|
|
|
|
const auto& gridView = simulator_.gridView();
|
|
|
|
const int number_of_cells = gridView.size(0);
|
|
|
|
|
|
|
|
cartesian_to_compressed_ = cartesianToCompressed(number_of_cells,
|
|
|
|
Opm::UgGridHelpers::globalCell(ugrid));
|
|
|
|
|
2018-10-16 08:59:16 -05:00
|
|
|
for (size_t i = 0; i < aquifersData.size(); ++i)
|
|
|
|
{
|
|
|
|
aquifers_Fetkovich.push_back(
|
2019-01-31 08:26:09 -06:00
|
|
|
AquiferFetkovich<TypeTag> (aquifer_connection.at(i), cartesian_to_compressed_, this->simulator_ , aquifersData.at(i))
|
2018-10-16 08:59:16 -05:00
|
|
|
);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<typename TypeTag>
|
|
|
|
bool
|
|
|
|
BlackoilAquiferModel<TypeTag>:: aquiferActive() const
|
|
|
|
{
|
|
|
|
return (aquiferCarterTracyActive() || aquiferFetkovichActive());
|
|
|
|
}
|
|
|
|
template<typename TypeTag>
|
|
|
|
bool
|
|
|
|
BlackoilAquiferModel<TypeTag>:: aquiferCarterTracyActive() const
|
|
|
|
{
|
|
|
|
return !aquifers_CarterTracy.empty();
|
|
|
|
}
|
|
|
|
template<typename TypeTag>
|
|
|
|
bool
|
|
|
|
BlackoilAquiferModel<TypeTag>:: aquiferFetkovichActive() const
|
|
|
|
{
|
|
|
|
return !aquifers_Fetkovich.empty();
|
|
|
|
}
|
|
|
|
} // namespace Opm
|