mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Get compressed to cartesian mapping and depths from the vanguard
This commit is contained in:
@@ -51,10 +51,9 @@ public:
|
||||
using Base::waterCompIdx;
|
||||
using Base::waterPhaseIdx;
|
||||
AquiferCarterTracy(const std::vector<Aquancon::AquancCell>& connections,
|
||||
const std::unordered_map<int, int>& cartesian_to_compressed,
|
||||
const Simulator& ebosSimulator,
|
||||
const AquiferCT::AQUCT_data& aquct_data)
|
||||
: Base(aquct_data.aquiferID, connections, cartesian_to_compressed, ebosSimulator)
|
||||
: Base(aquct_data.aquiferID, connections, ebosSimulator)
|
||||
, aquct_data_(aquct_data)
|
||||
{
|
||||
}
|
||||
@@ -92,10 +91,6 @@ protected:
|
||||
// This function is used to initialize and calculate the alpha_i for each grid connection to the aquifer
|
||||
inline void initializeConnections() override
|
||||
{
|
||||
const auto& eclState = this->ebos_simulator_.vanguard().eclState();
|
||||
const auto& grid = eclState.getInputGrid();
|
||||
|
||||
// We hack the cell depth values for now. We can actually get it from elementcontext pos
|
||||
this->cell_depth_.resize(this->size(), this->aquiferDepth());
|
||||
this->alphai_.resize(this->size(), 1.0);
|
||||
this->faceArea_connected_.resize(this->size(), 0.0);
|
||||
@@ -108,10 +103,13 @@ protected:
|
||||
this->cellToConnectionIdx_.resize(this->ebos_simulator_.gridView().size(/*codim=*/0), -1);
|
||||
for (size_t idx = 0; idx < this->size(); ++idx) {
|
||||
const auto global_index = this->connections_[idx].global_index;
|
||||
const int cell_index = this->cartesian_to_compressed_.at(global_index);
|
||||
const int cell_index = this->ebos_simulator_.vanguard().compressedIndex(global_index);
|
||||
|
||||
if (cell_index < 0) //the global_index is not part of this grid
|
||||
continue;
|
||||
|
||||
this->cellToConnectionIdx_[cell_index] = idx;
|
||||
const auto cellCenter = grid.getCellCenter(global_index);
|
||||
this->cell_depth_.at(idx) = cellCenter[2];
|
||||
this->cell_depth_.at(idx) = this->ebos_simulator_.vanguard().cellCenterDepth(cell_index);
|
||||
}
|
||||
// get default areas for all intersections
|
||||
const auto& gridView = this->ebos_simulator_.vanguard().gridView();
|
||||
|
||||
@@ -53,10 +53,9 @@ public:
|
||||
using Base::waterPhaseIdx;
|
||||
|
||||
AquiferFetkovich(const std::vector<Aquancon::AquancCell>& connections,
|
||||
const std::unordered_map<int, int>& cartesian_to_compressed,
|
||||
const Simulator& ebosSimulator,
|
||||
const Aquifetp::AQUFETP_data& aqufetp_data)
|
||||
: Base(aqufetp_data.aquiferID, connections, cartesian_to_compressed, ebosSimulator)
|
||||
: Base(aqufetp_data.aquiferID, connections, ebosSimulator)
|
||||
, aqufetp_data_(aqufetp_data)
|
||||
{
|
||||
}
|
||||
@@ -95,10 +94,6 @@ protected:
|
||||
|
||||
inline void initializeConnections() override
|
||||
{
|
||||
const auto& eclState = this->ebos_simulator_.vanguard().eclState();
|
||||
const auto& grid = eclState.getInputGrid();
|
||||
|
||||
// We hack the cell depth values for now. We can actually get it from elementcontext pos
|
||||
this->cell_depth_.resize(this->size(), this->aquiferDepth());
|
||||
this->alphai_.resize(this->size(), 1.0);
|
||||
this->faceArea_connected_.resize(this->size(), 0.0);
|
||||
@@ -111,10 +106,12 @@ protected:
|
||||
this->cellToConnectionIdx_.resize(this->ebos_simulator_.gridView().size(/*codim=*/0), -1);
|
||||
for (size_t idx = 0; idx < this->size(); ++idx) {
|
||||
const auto global_index = this->connections_[idx].global_index;
|
||||
const int cell_index = this->cartesian_to_compressed_.at(global_index);
|
||||
const int cell_index = this->ebos_simulator_.vanguard().compressedIndex(global_index);
|
||||
if (cell_index < 0) //the global_index is not part of this grid
|
||||
continue;
|
||||
|
||||
this->cellToConnectionIdx_[cell_index] = idx;
|
||||
const auto cellCenter = grid.getCellCenter(global_index);
|
||||
this->cell_depth_.at(idx) = cellCenter[2];
|
||||
this->cell_depth_.at(idx) = this->ebos_simulator_.vanguard().cellCenterDepth(cell_index);
|
||||
}
|
||||
// get areas for all connections
|
||||
const auto& gridView = this->ebos_simulator_.vanguard().gridView();
|
||||
|
||||
@@ -76,12 +76,10 @@ public:
|
||||
// Constructor
|
||||
AquiferInterface(int aqID,
|
||||
const std::vector<Aquancon::AquancCell>& connections,
|
||||
const std::unordered_map<int, int>& cartesian_to_compressed,
|
||||
const Simulator& ebosSimulator)
|
||||
: aquiferID(aqID)
|
||||
, connections_(connections)
|
||||
, ebos_simulator_(ebosSimulator)
|
||||
, cartesian_to_compressed_(cartesian_to_compressed)
|
||||
{
|
||||
}
|
||||
|
||||
@@ -216,7 +214,6 @@ protected:
|
||||
const int aquiferID;
|
||||
const std::vector<Aquancon::AquancCell> connections_;
|
||||
const Simulator& ebos_simulator_;
|
||||
const std::unordered_map<int, int> cartesian_to_compressed_;
|
||||
|
||||
// Grid variables
|
||||
std::vector<Scalar> faceArea_connected_;
|
||||
|
||||
@@ -83,7 +83,6 @@ protected:
|
||||
|
||||
Simulator& simulator_;
|
||||
|
||||
std::unordered_map<int, int> cartesian_to_compressed_;
|
||||
// TODO: probably we can use one variable to store both types of aquifers, because
|
||||
// they share the base class
|
||||
mutable std::vector<AquiferCarterTracy_object> aquifers_CarterTracy;
|
||||
|
||||
@@ -18,7 +18,6 @@
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <opm/grid/utility/cartesianToCompressed.hpp>
|
||||
namespace Opm
|
||||
{
|
||||
|
||||
@@ -170,20 +169,15 @@ BlackoilAquiferModel<TypeTag>::init()
|
||||
throw std::runtime_error("Aquifers currently do not work in parallel.");
|
||||
|
||||
// Get all the carter tracy aquifer properties data and put it in aquifers vector
|
||||
const int number_of_cells = simulator_.gridView().size(0);
|
||||
|
||||
cartesian_to_compressed_ = cartesianToCompressed(number_of_cells,
|
||||
this->simulator_.vanguard().grid().globalCell().data());
|
||||
|
||||
const auto& connections = aquifer.connections();
|
||||
for (const auto& aq : aquifer.ct()) {
|
||||
aquifers_CarterTracy.emplace_back(connections[aq.aquiferID],
|
||||
cartesian_to_compressed_, this->simulator_, aq);
|
||||
this->simulator_, aq);
|
||||
}
|
||||
|
||||
for (const auto& aq : aquifer.fetp()) {
|
||||
aquifers_Fetkovich.emplace_back(connections[aq.aquiferID],
|
||||
cartesian_to_compressed_, this->simulator_, aq);
|
||||
this->simulator_, aq);
|
||||
}
|
||||
}
|
||||
template <typename TypeTag>
|
||||
|
||||
Reference in New Issue
Block a user