Add Way of Introducing New Well Connections at Runtime

This commit adds a new member function

    void WellInterfaceGeneric<>::addPerforation()

which is, initially, intended as a way of introducing new well
connections arising from hydraulic fracturing.  We expect more work
in this area.
This commit is contained in:
Halvor M Nilsen 2025-01-31 18:31:56 +01:00 committed by Bård Skaflestad
parent 646d7c0a04
commit 9e235def1c
4 changed files with 82 additions and 1 deletions

View File

@ -1029,6 +1029,7 @@ list (APPEND PUBLIC_HEADER_FILES
opm/simulators/wells/RatioCalculator.hpp
opm/simulators/wells/RegionAttributeHelpers.hpp
opm/simulators/wells/RegionAverageCalculator.hpp
opm/simulators/wells/RuntimePerforation.hpp
opm/simulators/wells/SingleWellState.hpp
opm/simulators/wells/StandardWell.hpp
opm/simulators/wells/StandardWell_impl.hpp

View File

@ -0,0 +1,42 @@
/*
Copyright 2025 Equinor ASA.
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_RUNTIME_PERFORATION_HPP_INCLUDED
#define OPM_RUNTIME_PERFORATION_HPP_INCLUDED
namespace Opm {
/// Simple model of a well connection created at runtime, possibly as a
/// result of a geo-mechanical fracturing process.
struct RuntimePerforation
{
/// Active cell index, on current rank, that is dynamically perforated
/// by a well.
int cell{};
/// Connection's transmissibility factor.
double ctf{};
/// Depth at which the new connection is created.
double depth{};
};
} // namespace Opm
#endif // OPM_RUNTIME_PERFORATION_HPP_INCLUDED

View File

@ -20,6 +20,7 @@
*/
#include <config.h>
#include <opm/simulators/wells/WellInterfaceGeneric.hpp>
#include <opm/common/ErrorMacros.hpp>
@ -33,7 +34,9 @@
#include <opm/input/eclipse/Schedule/Well/WellPolymerProperties.hpp>
#include <opm/input/eclipse/Schedule/Well/WellTestState.hpp>
#include <opm/input/eclipse/Schedule/Well/WVFPEXP.hpp>
#include <opm/simulators/utils/DeferredLoggingErrorHelpers.hpp>
#include <opm/simulators/wells/PerforationData.hpp>
#include <opm/simulators/wells/ParallelWellInfo.hpp>
#include <opm/simulators/wells/VFPHelpers.hpp>
@ -45,10 +48,13 @@
#include <fmt/format.h>
#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstddef>
#include <stdexcept>
#include <string>
#include <vector>
namespace Opm {
@ -676,6 +682,36 @@ void WellInterfaceGeneric<Scalar>::resetWellOperability()
this->operability_status_.resetOperability();
}
template<class Scalar>
void WellInterfaceGeneric<Scalar>::addPerforations(const std::vector<RuntimePerforation>& perfs)
{
for (const auto& perf : perfs) {
auto it = std::find(well_cells_.begin(), well_cells_.end(), perf.cell);
if (it != this->well_cells_.end()) {
// If perforation to cell already exists, just add contribution.
const auto ind = std::distance(this->well_cells_.begin(), it);
this->well_index_[ind] += static_cast<Scalar>(perf.ctf);
}
else {
this->well_cells_.push_back(perf.cell);
this->well_index_.push_back(static_cast<Scalar>(perf.ctf));
this->perf_depth_.push_back(static_cast<Scalar>(perf.depth));
// Not strictly needed.
const double nan = std::nan("1");
this->perf_rep_radius_.push_back(nan);
this->perf_length_.push_back(nan);
this->bore_diameters_.push_back(nan);
// For now use the saturation table for the first cell.
this->saturation_table_number_
.push_back(this->saturation_table_number_.front());
++this->number_of_perforations_;
}
}
}
template<class Scalar>
Scalar WellInterfaceGeneric<Scalar>::wmicrobes_() const
{

View File

@ -26,6 +26,7 @@
#include <opm/input/eclipse/Schedule/Well/Well.hpp>
#include <opm/simulators/flow/BlackoilModelParameters.hpp>
#include <opm/simulators/wells/RuntimePerforation.hpp>
#include <map>
#include <optional>
@ -189,7 +190,6 @@ public:
void resetWellOperability();
virtual std::vector<Scalar> getPrimaryVars() const
{
return {};
@ -203,6 +203,8 @@ public:
virtual Scalar connectionDensity(const int globalConnIdx,
const int openConnIdx) const = 0;
void addPerforations(const std::vector<RuntimePerforation>& perfs);
protected:
bool getAllowCrossFlow() const;