mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
changed: move WBP BlackoilWellModel code to separate class
model using has-a instead of is-a in BlackoilWellModelGeneric
This commit is contained in:
parent
17d94b8054
commit
73ede837bb
@ -165,6 +165,7 @@ list (APPEND MAIN_SOURCE_FILES
|
||||
opm/simulators/wells/BlackoilWellModelGeneric.cpp
|
||||
opm/simulators/wells/BlackoilWellModelGuideRates.cpp
|
||||
opm/simulators/wells/BlackoilWellModelRestart.cpp
|
||||
opm/simulators/wells/BlackoilWellModelWBP.cpp
|
||||
opm/simulators/wells/ConnFiltrateData.cpp
|
||||
opm/simulators/wells/FractionCalculator.cpp
|
||||
opm/simulators/wells/GasLiftCommon.cpp
|
||||
@ -964,6 +965,8 @@ list (APPEND PUBLIC_HEADER_FILES
|
||||
opm/simulators/wells/BlackoilWellModelGeneric.hpp
|
||||
opm/simulators/wells/BlackoilWellModelGuideRates.hpp
|
||||
opm/simulators/wells/BlackoilWellModelRestart.hpp
|
||||
opm/simulators/wells/BlackoilWellModelWBP.hpp
|
||||
opm/simulators/wells/ConnectionIndexMap.hpp
|
||||
opm/simulators/wells/ConnFiltrateData.hpp
|
||||
opm/simulators/wells/FractionCalculator.hpp
|
||||
opm/simulators/wells/GasLiftCommon.hpp
|
||||
|
@ -237,7 +237,7 @@ template<class Scalar> class WellContributions;
|
||||
|
||||
data::WellBlockAveragePressures wellBlockAveragePressures() const
|
||||
{
|
||||
return this->computeWellBlockAveragePressures(this->gravity_);
|
||||
return this->wbp_.computeWellBlockAveragePressures(this->gravity_);
|
||||
}
|
||||
|
||||
#if COMPILE_GPU_BRIDGE
|
||||
|
@ -94,10 +94,10 @@ BlackoilWellModelGeneric(Schedule& schedule,
|
||||
, summaryState_(summaryState)
|
||||
, eclState_(eclState)
|
||||
, comm_(comm)
|
||||
, wbp_(*this)
|
||||
, phase_usage_(phase_usage)
|
||||
, terminal_output_(comm_.rank() == 0 &&
|
||||
Parameters::Get<Parameters::EnableTerminalOutput>())
|
||||
, wbpCalculationService_ { eclState.gridDims(), comm_ }
|
||||
, guideRate_(schedule)
|
||||
, active_wgstate_(phase_usage)
|
||||
, last_valid_wgstate_(phase_usage)
|
||||
@ -2033,138 +2033,6 @@ setPrimaryVarsDomain(const int domainIdx, const std::vector<Scalar>& vars)
|
||||
assert(offset == vars.size());
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelGeneric<Scalar>::
|
||||
registerOpenWellsForWBPCalculation()
|
||||
{
|
||||
assert (this->wbpCalcMap_.size() == this->wells_ecl_.size());
|
||||
|
||||
for (auto& wbpCalc : this->wbpCalcMap_) {
|
||||
wbpCalc.openWellIdx_.reset();
|
||||
}
|
||||
|
||||
auto openWellIdx = typename std::vector<WellInterfaceGeneric<Scalar>*>::size_type{0};
|
||||
for (const auto* openWell : this->well_container_generic_) {
|
||||
this->wbpCalcMap_[openWell->indexOfWell()].openWellIdx_ = openWellIdx++;
|
||||
}
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
typename ParallelWBPCalculation<Scalar>::EvaluatorFactory
|
||||
BlackoilWellModelGeneric<Scalar>::
|
||||
makeWellSourceEvaluatorFactory(const std::vector<Well>::size_type wellIdx) const
|
||||
{
|
||||
using Span = typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
|
||||
using Item = typename Span::Item;
|
||||
|
||||
return [wellIdx, this]() -> typename ParallelWBPCalculation<Scalar>::Evaluator
|
||||
{
|
||||
if (! this->wbpCalcMap_[wellIdx].openWellIdx_.has_value()) {
|
||||
// Well is stopped/shut. Return evaluator for stopped wells.
|
||||
return []([[maybe_unused]] const int connIdx, Span sourceTerm)
|
||||
{
|
||||
// Well/connection is stopped/shut. Set all items to
|
||||
// zero.
|
||||
|
||||
sourceTerm
|
||||
.set(Item::Pressure , 0.0)
|
||||
.set(Item::PoreVol , 0.0)
|
||||
.set(Item::MixtureDensity, 0.0)
|
||||
.set(Item::Depth , 0.0)
|
||||
;
|
||||
};
|
||||
}
|
||||
|
||||
// Well is open. Return an evaluator for open wells/open connections.
|
||||
return [this, wellPtr = this->well_container_generic_[*this->wbpCalcMap_[wellIdx].openWellIdx_]]
|
||||
(const int connIdx, Span sourceTerm)
|
||||
{
|
||||
// Note: The only item which actually matters for the WBP
|
||||
// calculation at the well reservoir connection level is the
|
||||
// mixture density. Set other items to zero.
|
||||
|
||||
const auto& connIdxMap =
|
||||
this->conn_idx_map_[wellPtr->indexOfWell()];
|
||||
|
||||
const auto rho = wellPtr->
|
||||
connectionDensity(connIdxMap.global(connIdx),
|
||||
connIdxMap.open(connIdx));
|
||||
|
||||
sourceTerm
|
||||
.set(Item::Pressure , 0.0)
|
||||
.set(Item::PoreVol , 0.0)
|
||||
.set(Item::MixtureDensity, rho)
|
||||
.set(Item::Depth , 0.0)
|
||||
;
|
||||
};
|
||||
};
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelGeneric<Scalar>::
|
||||
initializeWBPCalculationService()
|
||||
{
|
||||
this->wbpCalcMap_.clear();
|
||||
this->wbpCalcMap_.resize(this->wells_ecl_.size());
|
||||
|
||||
this->registerOpenWellsForWBPCalculation();
|
||||
|
||||
auto wellID = std::size_t{0};
|
||||
for (const auto& well : this->wells_ecl_) {
|
||||
this->wbpCalcMap_[wellID].wbpCalcIdx_ = this->wbpCalculationService_
|
||||
.createCalculator(well,
|
||||
this->local_parallel_well_info_[wellID],
|
||||
this->conn_idx_map_[wellID].local(),
|
||||
this->makeWellSourceEvaluatorFactory(wellID));
|
||||
|
||||
++wellID;
|
||||
}
|
||||
|
||||
this->wbpCalculationService_.defineCommunication();
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
data::WellBlockAveragePressures
|
||||
BlackoilWellModelGeneric<Scalar>::
|
||||
computeWellBlockAveragePressures(const Scalar gravity) const
|
||||
{
|
||||
auto wbpResult = data::WellBlockAveragePressures{};
|
||||
|
||||
using Calculated = typename PAvgCalculatorResult<Scalar>::WBPMode;
|
||||
using Output = data::WellBlockAvgPress::Quantity;
|
||||
|
||||
this->wbpCalculationService_.collectDynamicValues();
|
||||
|
||||
const auto numWells = this->wells_ecl_.size();
|
||||
for (auto wellID = 0*numWells; wellID < numWells; ++wellID) {
|
||||
const auto calcIdx = this->wbpCalcMap_[wellID].wbpCalcIdx_;
|
||||
const auto& well = this->wells_ecl_[wellID];
|
||||
|
||||
if (! well.hasRefDepth()) {
|
||||
// Can't perform depth correction without at least a
|
||||
// fall-back datum depth.
|
||||
continue;
|
||||
}
|
||||
|
||||
this->wbpCalculationService_
|
||||
.inferBlockAveragePressures(calcIdx, well.pavg(),
|
||||
gravity,
|
||||
well.getWPaveRefDepth());
|
||||
|
||||
const auto& result = this->wbpCalculationService_
|
||||
.averagePressures(calcIdx);
|
||||
|
||||
auto& reported = wbpResult.values[well.name()];
|
||||
|
||||
reported[Output::WBP] = result.value(Calculated::WBP);
|
||||
reported[Output::WBP4] = result.value(Calculated::WBP4);
|
||||
reported[Output::WBP5] = result.value(Calculated::WBP5);
|
||||
reported[Output::WBP9] = result.value(Calculated::WBP9);
|
||||
}
|
||||
|
||||
return wbpResult;
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelGeneric<Scalar>::
|
||||
reportGroupSwitching(DeferredLogger& local_deferredLogger) const
|
||||
|
@ -33,6 +33,8 @@
|
||||
|
||||
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
|
||||
|
||||
#include <opm/simulators/wells/BlackoilWellModelWBP.hpp>
|
||||
#include <opm/simulators/wells/ConnectionIndexMap.hpp>
|
||||
#include <opm/simulators/wells/ParallelPAvgDynamicSourceData.hpp>
|
||||
#include <opm/simulators/wells/ParallelWBPCalculation.hpp>
|
||||
#include <opm/simulators/wells/PerforationData.hpp>
|
||||
@ -111,6 +113,7 @@ public:
|
||||
// whether there exists any multisegment well open on this process
|
||||
bool anyMSWellOpenLocal() const;
|
||||
|
||||
const std::vector<Well>& eclWells() const { return wells_ecl_; }
|
||||
const Well& getWellEcl(const std::string& well_name) const;
|
||||
std::vector<Well> getLocalWells(const int timeStepIdx) const;
|
||||
const Schedule& schedule() const { return schedule_; }
|
||||
@ -267,6 +270,13 @@ public:
|
||||
this->closed_offending_wells_ == rhs.closed_offending_wells_;
|
||||
}
|
||||
|
||||
const ParallelWellInfo<Scalar>&
|
||||
parallelWellInfo(const std::size_t idx) const
|
||||
{ return local_parallel_well_info_[idx].get(); }
|
||||
|
||||
const ConnectionIndexMap& connectionIndexMap(const std::size_t idx)
|
||||
{ return conn_idx_map_[idx]; }
|
||||
|
||||
protected:
|
||||
/*
|
||||
The dynamic state of the well model is maintained with an instance
|
||||
@ -453,20 +463,11 @@ protected:
|
||||
void assignMassGasRate(data::Wells& wsrpt,
|
||||
const Scalar& gasDensity) const;
|
||||
|
||||
void registerOpenWellsForWBPCalculation();
|
||||
|
||||
typename ParallelWBPCalculation<Scalar>::EvaluatorFactory
|
||||
makeWellSourceEvaluatorFactory(const std::vector<Well>::size_type wellIdx) const;
|
||||
|
||||
void initializeWBPCalculationService();
|
||||
|
||||
data::WellBlockAveragePressures
|
||||
computeWellBlockAveragePressures(const Scalar gravity) const;
|
||||
|
||||
Schedule& schedule_;
|
||||
const SummaryState& summaryState_;
|
||||
const EclipseState& eclState_;
|
||||
const Parallel::Communication& comm_;
|
||||
BlackoilWellModelWBP<Scalar> wbp_;
|
||||
|
||||
PhaseUsage phase_usage_;
|
||||
bool terminal_output_{false};
|
||||
@ -486,92 +487,6 @@ protected:
|
||||
// Times at which wells were shut (for WCYCLE)
|
||||
std::map<std::string, double> well_close_times_;
|
||||
|
||||
/// Connection index mappings
|
||||
class ConnectionIndexMap
|
||||
{
|
||||
public:
|
||||
/// Constructor.
|
||||
///
|
||||
/// \param[in] numConns Total number of well connections, both open
|
||||
/// and closed/shut. Typically \code WellConnections::size() \endcode.
|
||||
explicit ConnectionIndexMap(const std::size_t numConns)
|
||||
: local_(numConns, -1)
|
||||
{
|
||||
this->global_.reserve(numConns);
|
||||
this->open_.reserve(numConns);
|
||||
}
|
||||
|
||||
/// Enumerate/map new active connection.
|
||||
///
|
||||
/// \param[in] connIdx Global well connection index. Must be an
|
||||
/// integer in the range 0..numConns-1.
|
||||
///
|
||||
/// \param[in] connIsOpen Whether or not the connection is
|
||||
/// open/flowing.
|
||||
void addActiveConnection(const int connIdx,
|
||||
const bool connIsOpen)
|
||||
{
|
||||
this->local_[connIdx] =
|
||||
static_cast<int>(this->global_.size());
|
||||
|
||||
this->global_.push_back(connIdx);
|
||||
|
||||
const auto open_conn_idx = connIsOpen
|
||||
? this->num_open_conns_++
|
||||
: -1;
|
||||
|
||||
this->open_.push_back(open_conn_idx);
|
||||
}
|
||||
|
||||
/// Get local connection IDs/indices of every existing well
|
||||
/// connection.
|
||||
///
|
||||
/// Negative value (-1) for connections that don't intersect the
|
||||
/// current rank.
|
||||
const std::vector<int>& local() const
|
||||
{
|
||||
return this->local_;
|
||||
}
|
||||
|
||||
/// Get global connection ID of local (on-rank) connection.
|
||||
///
|
||||
/// \param[in] connIdx Local connection index.
|
||||
///
|
||||
/// \return Global connection ID of \p connIdx.
|
||||
int global(const int connIdx) const
|
||||
{
|
||||
return this->global_[connIdx];
|
||||
}
|
||||
|
||||
/// Get open connection ID of local (on-rank) connection.
|
||||
///
|
||||
/// \param[in] connIdx Local connection index.
|
||||
///
|
||||
/// \return Open connection ID of \p connIdx. Integer in the range
|
||||
/// 0..#open connections - 1 if the connection is open or negative
|
||||
/// value (-1) otherwise.
|
||||
int open(const int connIdx) const
|
||||
{
|
||||
return this->open_[connIdx];
|
||||
}
|
||||
|
||||
private:
|
||||
/// Local connection IDs/indices of every existing well connection.
|
||||
/// Negative value (-1) for connections that don't intersect the
|
||||
/// current rank.
|
||||
std::vector<int> local_{};
|
||||
|
||||
/// Global connection index of each on-rank reservoir connection.
|
||||
/// Reverse/transpose mapping of \c local_.
|
||||
std::vector<int> global_{};
|
||||
|
||||
/// Open connection index of each on-rank reservoir connection.
|
||||
std::vector<int> open_{};
|
||||
|
||||
/// Number of open connections on this rank.
|
||||
int num_open_conns_{0};
|
||||
};
|
||||
|
||||
std::vector<ConnectionIndexMap> conn_idx_map_{};
|
||||
std::function<bool(const Well&)> not_on_process_{};
|
||||
|
||||
@ -584,15 +499,6 @@ protected:
|
||||
std::vector<std::reference_wrapper<ParallelWellInfo<Scalar>>> local_parallel_well_info_;
|
||||
|
||||
std::vector<WellProdIndexCalculator<Scalar>> prod_index_calc_;
|
||||
mutable ParallelWBPCalculation<Scalar> wbpCalculationService_;
|
||||
|
||||
struct WBPCalcID
|
||||
{
|
||||
std::optional<typename std::vector<WellInterfaceGeneric<Scalar>*>::size_type> openWellIdx_{};
|
||||
std::size_t wbpCalcIdx_{};
|
||||
};
|
||||
|
||||
std::vector<WBPCalcID> wbpCalcMap_{};
|
||||
|
||||
std::vector<int> pvt_region_idx_;
|
||||
|
||||
|
188
opm/simulators/wells/BlackoilWellModelWBP.cpp
Normal file
188
opm/simulators/wells/BlackoilWellModelWBP.cpp
Normal file
@ -0,0 +1,188 @@
|
||||
/*
|
||||
Copyright 2016 SINTEF ICT, Applied Mathematics.
|
||||
Copyright 2016 - 2017 Statoil ASA.
|
||||
Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
|
||||
Copyright 2016 - 2018 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#include <config.h>
|
||||
#include <opm/simulators/wells/BlackoilWellModelWBP.hpp>
|
||||
|
||||
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
|
||||
|
||||
#include <opm/simulators/wells/BlackoilWellModelGeneric.hpp>
|
||||
|
||||
#include <cassert>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
template<class Scalar>
|
||||
BlackoilWellModelWBP<Scalar>::
|
||||
BlackoilWellModelWBP(BlackoilWellModelGeneric<Scalar>& well_model)
|
||||
: well_model_(well_model)
|
||||
, wbpCalculationService_(well_model.eclipseState().gridDims(), well_model.comm())
|
||||
{}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelWBP<Scalar>::
|
||||
initializeSources(typename ParallelWBPCalculation<Scalar>::GlobalToLocal index,
|
||||
typename ParallelWBPCalculation<Scalar>::Evaluator eval)
|
||||
{
|
||||
this->wbpCalculationService_.localCellIndex(index).evalCellSource(eval);
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelWBP<Scalar>::
|
||||
registerOpenWellsForWBPCalculation()
|
||||
{
|
||||
assert(this->wbpCalcMap_.size() ==
|
||||
static_cast<std::size_t>(well_model_.numLocalWells()));
|
||||
|
||||
for (auto& wbpCalc : this->wbpCalcMap_) {
|
||||
wbpCalc.openWellIdx_.reset();
|
||||
}
|
||||
|
||||
auto openWellIdx = typename std::vector<WellInterfaceGeneric<Scalar>*>::size_type{0};
|
||||
for (const auto* openWell : well_model_.genericWells()) {
|
||||
this->wbpCalcMap_[openWell->indexOfWell()].openWellIdx_ = openWellIdx++;
|
||||
}
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
void BlackoilWellModelWBP<Scalar>::
|
||||
initializeWBPCalculationService()
|
||||
{
|
||||
this->wbpCalcMap_.clear();
|
||||
this->wbpCalcMap_.resize(well_model_.numLocalWells());
|
||||
|
||||
this->registerOpenWellsForWBPCalculation();
|
||||
|
||||
auto wellID = std::size_t{0};
|
||||
for (const auto& well : well_model_.eclWells()) {
|
||||
this->wbpCalcMap_[wellID].wbpCalcIdx_ = this->wbpCalculationService_
|
||||
.createCalculator(well,
|
||||
well_model_.parallelWellInfo(wellID),
|
||||
well_model_.connectionIndexMap(wellID).local(),
|
||||
this->makeWellSourceEvaluatorFactory(wellID));
|
||||
|
||||
++wellID;
|
||||
}
|
||||
|
||||
this->wbpCalculationService_.defineCommunication();
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
data::WellBlockAveragePressures
|
||||
BlackoilWellModelWBP<Scalar>::
|
||||
computeWellBlockAveragePressures(const Scalar gravity) const
|
||||
{
|
||||
auto wbpResult = data::WellBlockAveragePressures{};
|
||||
|
||||
using Calculated = typename PAvgCalculatorResult<Scalar>::WBPMode;
|
||||
using Output = data::WellBlockAvgPress::Quantity;
|
||||
|
||||
this->wbpCalculationService_.collectDynamicValues();
|
||||
|
||||
const auto numWells = well_model_.numLocalWells();
|
||||
for (auto wellID = 0*numWells; wellID < numWells; ++wellID) {
|
||||
const auto calcIdx = this->wbpCalcMap_[wellID].wbpCalcIdx_;
|
||||
const auto& well = well_model_.eclWells()[wellID];
|
||||
|
||||
if (! well.hasRefDepth()) {
|
||||
// Can't perform depth correction without at least a
|
||||
// fall-back datum depth.
|
||||
continue;
|
||||
}
|
||||
|
||||
this->wbpCalculationService_
|
||||
.inferBlockAveragePressures(calcIdx, well.pavg(),
|
||||
gravity,
|
||||
well.getWPaveRefDepth());
|
||||
|
||||
const auto& result = this->wbpCalculationService_
|
||||
.averagePressures(calcIdx);
|
||||
|
||||
auto& reported = wbpResult.values[well.name()];
|
||||
|
||||
reported[Output::WBP] = result.value(Calculated::WBP);
|
||||
reported[Output::WBP4] = result.value(Calculated::WBP4);
|
||||
reported[Output::WBP5] = result.value(Calculated::WBP5);
|
||||
reported[Output::WBP9] = result.value(Calculated::WBP9);
|
||||
}
|
||||
|
||||
return wbpResult;
|
||||
}
|
||||
|
||||
template<class Scalar>
|
||||
typename ParallelWBPCalculation<Scalar>::EvaluatorFactory
|
||||
BlackoilWellModelWBP<Scalar>::
|
||||
makeWellSourceEvaluatorFactory(const std::vector<Well>::size_type wellIdx) const
|
||||
{
|
||||
using Span = typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
|
||||
using Item = typename Span::Item;
|
||||
|
||||
return [wellIdx, this]() -> typename ParallelWBPCalculation<Scalar>::Evaluator
|
||||
{
|
||||
if (! this->wbpCalcMap_[wellIdx].openWellIdx_.has_value()) {
|
||||
// Well is stopped/shut. Return evaluator for stopped wells.
|
||||
return []([[maybe_unused]] const int connIdx, Span sourceTerm)
|
||||
{
|
||||
// Well/connection is stopped/shut. Set all items to
|
||||
// zero.
|
||||
|
||||
sourceTerm
|
||||
.set(Item::Pressure , 0.0)
|
||||
.set(Item::PoreVol , 0.0)
|
||||
.set(Item::MixtureDensity, 0.0)
|
||||
.set(Item::Depth , 0.0)
|
||||
;
|
||||
};
|
||||
}
|
||||
|
||||
// Well is open. Return an evaluator for open wells/open connections.
|
||||
return [this, wellPtr = well_model_.genericWells()[*this->wbpCalcMap_[wellIdx].openWellIdx_]]
|
||||
(const int connIdx, Span sourceTerm)
|
||||
{
|
||||
// Note: The only item which actually matters for the WBP
|
||||
// calculation at the well reservoir connection level is the
|
||||
// mixture density. Set other items to zero.
|
||||
|
||||
const auto& connIdxMap =
|
||||
well_model_.connectionIndexMap(wellPtr->indexOfWell());
|
||||
|
||||
const auto rho = wellPtr->
|
||||
connectionDensity(connIdxMap.global(connIdx),
|
||||
connIdxMap.open(connIdx));
|
||||
|
||||
sourceTerm
|
||||
.set(Item::Pressure , 0.0)
|
||||
.set(Item::PoreVol , 0.0)
|
||||
.set(Item::MixtureDensity, rho)
|
||||
.set(Item::Depth , 0.0)
|
||||
;
|
||||
};
|
||||
};
|
||||
}
|
||||
|
||||
template class BlackoilWellModelWBP<double>;
|
||||
|
||||
#if FLOW_INSTANTIATE_FLOAT
|
||||
template class BlackoilWellModelWBP<float>;
|
||||
#endif
|
||||
|
||||
}
|
75
opm/simulators/wells/BlackoilWellModelWBP.hpp
Normal file
75
opm/simulators/wells/BlackoilWellModelWBP.hpp
Normal file
@ -0,0 +1,75 @@
|
||||
/*
|
||||
Copyright 2016 SINTEF ICT, Applied Mathematics.
|
||||
Copyright 2016 - 2017 Statoil ASA.
|
||||
Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
|
||||
Copyright 2016 - 2018 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#ifndef OPM_BLACKOILWELLMODEL_WBP_HEADER_INCLUDED
|
||||
#define OPM_BLACKOILWELLMODEL_WBP_HEADER_INCLUDED
|
||||
|
||||
#include <opm/output/data/Wells.hpp>
|
||||
|
||||
#include <opm/simulators/wells/ParallelWBPCalculation.hpp>
|
||||
#include <opm/simulators/wells/WellInterfaceGeneric.hpp>
|
||||
|
||||
#include <cstddef>
|
||||
#include <optional>
|
||||
#include <vector>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
template<class Scalar> class BlackoilWellModelGeneric;
|
||||
|
||||
/// Class for handling the blackoil well model.
|
||||
template<class Scalar>
|
||||
class BlackoilWellModelWBP
|
||||
{
|
||||
public:
|
||||
BlackoilWellModelWBP(BlackoilWellModelGeneric<Scalar>& well_model);
|
||||
|
||||
void initializeSources(typename ParallelWBPCalculation<Scalar>::GlobalToLocal index,
|
||||
typename ParallelWBPCalculation<Scalar>::Evaluator eval);
|
||||
|
||||
void registerOpenWellsForWBPCalculation();
|
||||
|
||||
typename ParallelWBPCalculation<Scalar>::EvaluatorFactory
|
||||
makeWellSourceEvaluatorFactory(const std::vector<Well>::size_type wellIdx) const;
|
||||
|
||||
void initializeWBPCalculationService();
|
||||
|
||||
data::WellBlockAveragePressures
|
||||
computeWellBlockAveragePressures(const Scalar gravity) const;
|
||||
|
||||
private:
|
||||
BlackoilWellModelGeneric<Scalar>& well_model_;
|
||||
mutable ParallelWBPCalculation<Scalar> wbpCalculationService_;
|
||||
|
||||
struct WBPCalcID
|
||||
{
|
||||
std::optional<typename std::vector<WellInterfaceGeneric<Scalar>*>::size_type> openWellIdx_{};
|
||||
std::size_t wbpCalcIdx_{};
|
||||
};
|
||||
|
||||
std::vector<WBPCalcID> wbpCalcMap_{};
|
||||
};
|
||||
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif
|
@ -98,10 +98,10 @@ namespace Opm {
|
||||
using SourceDataSpan =
|
||||
typename PAvgDynamicSourceData<Scalar>::template SourceDataSpan<Scalar>;
|
||||
|
||||
this->wbpCalculationService_
|
||||
.localCellIndex([this](const std::size_t globalIndex)
|
||||
{ return this->compressedIndexForInterior(globalIndex); })
|
||||
.evalCellSource([this](const int localCell, SourceDataSpan sourceTerms)
|
||||
this->wbp_.initializeSources(
|
||||
[this](const std::size_t globalIndex)
|
||||
{ return this->compressedIndexForInterior(globalIndex); },
|
||||
[this](const int localCell, SourceDataSpan sourceTerms)
|
||||
{
|
||||
using Item = typename SourceDataSpan::Item;
|
||||
|
||||
@ -139,7 +139,8 @@ namespace Opm {
|
||||
if (haveWat) { rho += weightedPhaseDensity(iw); }
|
||||
|
||||
sourceTerms.set(Item::MixtureDensity, rho);
|
||||
});
|
||||
}
|
||||
);
|
||||
}
|
||||
|
||||
template<typename TypeTag>
|
||||
@ -261,7 +262,7 @@ namespace Opm {
|
||||
{
|
||||
this->initializeWellPerfData();
|
||||
this->initializeWellState(reportStepIdx);
|
||||
this->initializeWBPCalculationService();
|
||||
this->wbp_.initializeWBPCalculationService();
|
||||
|
||||
if (this->param_.use_multisegment_well_ && this->anyMSWellOpenLocal()) {
|
||||
this->wellState().initWellStateMSWell(this->wells_ecl_, &this->prevWellState());
|
||||
@ -961,7 +962,7 @@ namespace Opm {
|
||||
}
|
||||
}
|
||||
|
||||
this->registerOpenWellsForWBPCalculation();
|
||||
this->wbp_.registerOpenWellsForWBPCalculation();
|
||||
}
|
||||
|
||||
|
||||
|
119
opm/simulators/wells/ConnectionIndexMap.hpp
Normal file
119
opm/simulators/wells/ConnectionIndexMap.hpp
Normal file
@ -0,0 +1,119 @@
|
||||
/*
|
||||
Copyright 2016 SINTEF ICT, Applied Mathematics.
|
||||
Copyright 2016 - 2017 Statoil ASA.
|
||||
Copyright 2017 Dr. Blatt - HPC-Simulation-Software & Services
|
||||
Copyright 2016 - 2018 IRIS AS
|
||||
|
||||
This file is part of the Open Porous Media project (OPM).
|
||||
|
||||
OPM is free software: you can redistribute it and/or modify
|
||||
it under the terms of the GNU General Public License as published by
|
||||
the Free Software Foundation, either version 3 of the License, or
|
||||
(at your option) any later version.
|
||||
|
||||
OPM is distributed in the hope that it will be useful,
|
||||
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
||||
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
||||
GNU General Public License for more details.
|
||||
|
||||
You should have received a copy of the GNU General Public License
|
||||
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
||||
*/
|
||||
|
||||
#ifndef OPM_CONNECTION_INDEX_MAP_HEADER_INCLUDED
|
||||
#define OPM_CONNECTION_INDEX_MAP_HEADER_INCLUDED
|
||||
|
||||
#include <cstddef>
|
||||
#include <vector>
|
||||
|
||||
namespace Opm {
|
||||
|
||||
/// Connection index mappings
|
||||
class ConnectionIndexMap
|
||||
{
|
||||
public:
|
||||
/// Constructor.
|
||||
///
|
||||
/// \param[in] numConns Total number of well connections, both open
|
||||
/// and closed/shut. Typically \code WellConnections::size() \endcode.
|
||||
explicit ConnectionIndexMap(const std::size_t numConns)
|
||||
: local_(numConns, -1)
|
||||
{
|
||||
this->global_.reserve(numConns);
|
||||
this->open_.reserve(numConns);
|
||||
}
|
||||
|
||||
/// Enumerate/map new active connection.
|
||||
///
|
||||
/// \param[in] connIdx Global well connection index. Must be an
|
||||
/// integer in the range 0..numConns-1.
|
||||
///
|
||||
/// \param[in] connIsOpen Whether or not the connection is
|
||||
/// open/flowing.
|
||||
void addActiveConnection(const int connIdx,
|
||||
const bool connIsOpen)
|
||||
{
|
||||
this->local_[connIdx] =
|
||||
static_cast<int>(this->global_.size());
|
||||
|
||||
this->global_.push_back(connIdx);
|
||||
|
||||
const auto open_conn_idx = connIsOpen
|
||||
? this->num_open_conns_++
|
||||
: -1;
|
||||
|
||||
this->open_.push_back(open_conn_idx);
|
||||
}
|
||||
|
||||
/// Get local connection IDs/indices of every existing well
|
||||
/// connection.
|
||||
///
|
||||
/// Negative value (-1) for connections that don't intersect the
|
||||
/// current rank.
|
||||
const std::vector<int>& local() const
|
||||
{
|
||||
return this->local_;
|
||||
}
|
||||
|
||||
/// Get global connection ID of local (on-rank) connection.
|
||||
///
|
||||
/// \param[in] connIdx Local connection index.
|
||||
///
|
||||
/// \return Global connection ID of \p connIdx.
|
||||
int global(const int connIdx) const
|
||||
{
|
||||
return this->global_[connIdx];
|
||||
}
|
||||
|
||||
/// Get open connection ID of local (on-rank) connection.
|
||||
///
|
||||
/// \param[in] connIdx Local connection index.
|
||||
///
|
||||
/// \return Open connection ID of \p connIdx. Integer in the range
|
||||
/// 0..#open connections - 1 if the connection is open or negative
|
||||
/// value (-1) otherwise.
|
||||
int open(const int connIdx) const
|
||||
{
|
||||
return this->open_[connIdx];
|
||||
}
|
||||
|
||||
private:
|
||||
/// Local connection IDs/indices of every existing well connection.
|
||||
/// Negative value (-1) for connections that don't intersect the
|
||||
/// current rank.
|
||||
std::vector<int> local_{};
|
||||
|
||||
/// Global connection index of each on-rank reservoir connection.
|
||||
/// Reverse/transpose mapping of \c local_.
|
||||
std::vector<int> global_{};
|
||||
|
||||
/// Open connection index of each on-rank reservoir connection.
|
||||
std::vector<int> open_{};
|
||||
|
||||
/// Number of open connections on this rank.
|
||||
int num_open_conns_{0};
|
||||
};
|
||||
|
||||
} // namespace Opm
|
||||
|
||||
#endif
|
Loading…
Reference in New Issue
Block a user