Output rates to multisegment wells

This commit is contained in:
Svenn Tveit 2024-06-05 15:13:05 +02:00
parent db970d58d2
commit adc36d64a5
5 changed files with 52 additions and 2 deletions

View File

@ -89,6 +89,8 @@ public:
getWellFreeTracerRates() const {return wellFreeTracerRate_;}
const std::map<std::pair<std::string, std::string>, Scalar>&
getWellSolTracerRates() const {return wellSolTracerRate_;}
const std::map<std::tuple<std::string, std::string, std::size_t>, Scalar>&
getMswTracerRates() const {return mSwTracerRate_;}
template<class Serializer>
void serializeOp(Serializer& serializer)
@ -134,10 +136,14 @@ protected:
std::vector<TracerVectorSingle> freeTracerConcentration_;
std::vector<TracerVectorSingle> solTracerConcentration_;
// <wellName, tracerIdx> -> wellRate
// <wellName, tracerName> -> wellRate
std::map<std::pair<std::string, std::string>, Scalar> wellTracerRate_;
std::map<std::pair<std::string, std::string>, Scalar> wellFreeTracerRate_;
std::map<std::pair<std::string, std::string>, Scalar> wellSolTracerRate_;
// <wellName, tracerName, segNum> -> wellRate
std::map<std::tuple<std::string, std::string, std::size_t>, Scalar> mSwTracerRate_;
/// \brief Function returning the cell centers
std::function<std::array<double,dimWorld>(int)> centroids_;
};

View File

@ -489,7 +489,6 @@ protected:
}
Scalar rate_f = rate - rate_s;
if (rate_f > 0) {
for (int tIdx = 0; tIdx < tr.numTracer(); ++tIdx) {
// Injection of free tracer only
@ -498,6 +497,10 @@ protected:
// Store _injector_ tracer rate for reporting
this->wellTracerRate_.at(std::make_pair(eclWell.name(),this->name(tr.idx_[tIdx]))) += rate_f*wtracer[tIdx];
this->wellFreeTracerRate_.at(std::make_pair(eclWell.name(),this->wellfname(tr.idx_[tIdx]))) += rate_f*wtracer[tIdx];
if (eclWell.isMultiSegment()) {
this->mSwTracerRate_[std::make_tuple(eclWell.name(), this->name(tr.idx_[tIdx]), eclWell.getConnections().get(i).segment())] =
rate_f*wtracer[tIdx];
}
}
dfVol_[tr.phaseIdx_][I] -= rate_f * dt;
}
@ -511,6 +514,10 @@ protected:
rate_f * tr.concentration_[tIdx][I][0];
this->wellFreeTracerRate_.at(std::make_pair(eclWell.name(),this->wellfname(tr.idx_[tIdx]))) +=
rate_f * tr.concentration_[tIdx][I][0];
if (eclWell.isMultiSegment()) {
this->mSwTracerRate_[std::make_tuple(eclWell.name(), this->name(tr.idx_[tIdx]), eclWell.getConnections().get(i).segment())] =
rate_f * tr.concentration_[tIdx][I][0];
}
}
dfVol_[tr.phaseIdx_][I] -= rate_f * dt;
@ -525,6 +532,10 @@ protected:
rate_s * tr.concentration_[tIdx][I][1];
this->wellSolTracerRate_.at(std::make_pair(eclWell.name(),this->wellsname(tr.idx_[tIdx]))) +=
rate_s * tr.concentration_[tIdx][I][1];
if (eclWell.isMultiSegment()) {
this->mSwTracerRate_[std::make_tuple(eclWell.name(), this->name(tr.idx_[tIdx]), eclWell.getConnections().get(i).segment())] =
rate_s * tr.concentration_[tIdx][I][1];
}
}
dsVol_[tr.phaseIdx_][I] -= rate_s * dt;

View File

@ -271,10 +271,13 @@ template<class Scalar> class WellContributions;
.tracerModel().getWellFreeTracerRates();
const auto& solTracerRates = simulator_.problem()
.tracerModel().getWellSolTracerRates();
const auto& mswTracerRates = simulator_.problem()
.tracerModel().getMswTracerRates();
this->assignWellTracerRates(wsrpt, tracerRates);
this->assignWellTracerRates(wsrpt, freeTracerRates);
this->assignWellTracerRates(wsrpt, solTracerRates);
this->assignMswTracerRates(wsrpt, mswTracerRates);
}
BlackoilWellModelGuideRates(*this)

View File

@ -1698,6 +1698,33 @@ assignWellTracerRates(data::Wells& wsrpt,
}
}
template<class Scalar>
void BlackoilWellModelGeneric<Scalar>::
assignMswTracerRates(data::Wells& wsrpt,
const MswTracerRates& mswTracerRates) const
{
if (mswTracerRates.empty())
return;
for (const auto& mswTR : mswTracerRates) {
std::string wellName = std::get<0>(mswTR.first);
auto xwPos = wsrpt.find(wellName);
if (xwPos == wsrpt.end()) { // No well results.
continue;
}
std::string tracerName = std::get<1>(mswTR.first);
std::size_t segNumber = std::get<2>(mswTR.first);
Scalar rate = mswTR.second;
auto& wData = xwPos->second;
auto segPos = wData.segments.find(segNumber);
if (segPos != wData.segments.end()) {
auto& segment = segPos->second;
segment.rates.set(data::Rates::opt::tracer, rate, tracerName);
}
}
}
template<class Scalar>
std::vector<std::vector<int>> BlackoilWellModelGeneric<Scalar>::
getMaxWellConnections() const

View File

@ -437,6 +437,9 @@ protected:
using WellTracerRates = std::map<std::pair<std::string, std::string>, Scalar>;
void assignWellTracerRates(data::Wells& wsrpt,
const WellTracerRates& wellTracerRates) const;
using MswTracerRates = std::map<std::tuple<std::string, std::string, std::size_t>, Scalar>;
void assignMswTracerRates(data::Wells& wsrpt,
const MswTracerRates& mswTracerRates) const;
Schedule& schedule_;
const SummaryState& summaryState_;