mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
refining the function updateInjFCMult
for code improvements.
This commit is contained in:
parent
22269c92c3
commit
e264042c62
@ -315,7 +315,7 @@ namespace Opm {
|
||||
const auto it = this->filtration_particle_volume_.find(well->name());
|
||||
if (it != this->filtration_particle_volume_.end()) {
|
||||
const auto& filtration_particle_volume = it->second;
|
||||
well->updateInjFCMult(filtration_particle_volume);
|
||||
well->updateInjFCMult(filtration_particle_volume, local_deferredLogger);
|
||||
}
|
||||
}
|
||||
}
|
||||
|
@ -749,9 +749,7 @@ updateFiltrationParticleVolume(const double dt, const size_t water_index,
|
||||
}
|
||||
|
||||
void WellInterfaceGeneric::
|
||||
updateInjFCMult(const std::vector<double>& filtration_particle_volume) {
|
||||
// the filter injecting concentration, the porosity, the area size
|
||||
// we also need the permeability of the formation, and rw and some other information
|
||||
updateInjFCMult(const std::vector<double>& filtration_particle_volume, DeferredLogger& deferred_logger) {
|
||||
for (int perf = 0; perf < this->number_of_perforations_; ++perf) {
|
||||
const auto perf_ecl_index = this->perforationData()[perf].ecl_index;
|
||||
const auto& connections = this->well_ecl_.getConnections();
|
||||
@ -760,27 +758,38 @@ updateInjFCMult(const std::vector<double>& filtration_particle_volume) {
|
||||
const auto& filter_cake = connection.getFilterCake();
|
||||
const double area = connection.getFilterCakeArea();
|
||||
const double poro = filter_cake.poro;
|
||||
const double thickness = filtration_particle_volume[perf] / (area * (1. - poro));
|
||||
const double perm = filter_cake.perm;
|
||||
const double rw = connection.getFilterCakeRadius();
|
||||
const auto cr0 = connection.r0();
|
||||
const auto crw = connection.rw();
|
||||
const auto cskinfactor = connection.skinFactor();
|
||||
const double K = connection.Kh() / connection.connectionLength();
|
||||
const double factor = filter_cake.sf_multiplier;
|
||||
// compute a multiplier for the well connection transmissibility
|
||||
if (filter_cake.geometry == FilterCake::FilterCakeGeometry::LINEAR) {
|
||||
const double skin_factor = thickness / rw * K / perm * factor;
|
||||
const auto denom = std::log(cr0 / std::min(crw, cr0)) + cskinfactor;
|
||||
const auto denom2 = denom + skin_factor;
|
||||
this->inj_fc_multiplier_[perf] = denom / denom2;
|
||||
} else if (filter_cake.geometry == FilterCake::FilterCakeGeometry::RADIAL) {
|
||||
const double rc = std::sqrt(rw * rw + 2. * rw * thickness );
|
||||
const double skin_factor = K / perm * std::log(rc/rw) * factor;
|
||||
const auto denom = std::log(cr0 / std::min(crw, cr0)) + cskinfactor;
|
||||
const auto denom2 = denom + skin_factor;
|
||||
this->inj_fc_multiplier_[perf] = denom / denom2;
|
||||
// the thickness of the filtration cake
|
||||
const double thickness = filtration_particle_volume[perf] / (area * (1. - poro));
|
||||
|
||||
double skin_factor = 0.;
|
||||
switch (filter_cake.geometry) {
|
||||
case FilterCake::FilterCakeGeometry::LINEAR: {
|
||||
skin_factor = thickness / rw * K / perm * factor;
|
||||
break;
|
||||
}
|
||||
case FilterCake::FilterCakeGeometry::RADIAL: {
|
||||
const double rc = std::sqrt(rw * rw + 2. * rw * thickness);
|
||||
skin_factor = K / perm * std::log(rc / rw) * factor;
|
||||
break;
|
||||
}
|
||||
default:
|
||||
OPM_DEFLOG_THROW(std::runtime_error,
|
||||
fmt::format(" Invalid filtration cake geometry type ({}) for well {}",
|
||||
FilterCake::filterCakeGeometryToString(filter_cake.geometry), name()),
|
||||
deferred_logger);
|
||||
}
|
||||
// the original skin factor for the connection
|
||||
const auto cskinfactor = connection.skinFactor();
|
||||
// compute a multiplier for the well connection transmissibility
|
||||
const auto denom = std::log(cr0 / std::min(crw, cr0)) + cskinfactor;
|
||||
const auto denom2 = denom + skin_factor;
|
||||
this->inj_fc_multiplier_[perf] = denom / denom2;
|
||||
} else {
|
||||
this->inj_fc_multiplier_[perf] = 1.0;
|
||||
}
|
||||
|
@ -193,7 +193,7 @@ public:
|
||||
const WellState& well_state, std::vector<double>& filtration_particle_volume) const;
|
||||
|
||||
// update the multiplier for well transmissbility due to cake filteration
|
||||
void updateInjFCMult(const std::vector<double>& filtration_particle_volume);
|
||||
void updateInjFCMult(const std::vector<double>& filtration_particle_volume, DeferredLogger& deferred_logger);
|
||||
|
||||
// whether a well is specified with a non-zero and valid VFP table number
|
||||
bool isVFPActive(DeferredLogger& deferred_logger) const;
|
||||
|
Loading…
Reference in New Issue
Block a user