Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

changed: move NLDD BlackoilWellModel code to separate class #5847

Open
wants to merge 1 commit into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 3 additions & 0 deletions CMakeLists_files.cmake
Original file line number Diff line number Diff line change
Expand Up @@ -167,6 +167,7 @@ list (APPEND MAIN_SOURCE_FILES
opm/simulators/wells/BlackoilWellModelGasLift.cpp
opm/simulators/wells/BlackoilWellModelGeneric.cpp
opm/simulators/wells/BlackoilWellModelGuideRates.cpp
opm/simulators/wells/BlackoilWellModelNldd.cpp
opm/simulators/wells/BlackoilWellModelRestart.cpp
opm/simulators/wells/BlackoilWellModelWBP.cpp
opm/simulators/wells/ConnFiltrateData.cpp
Expand Down Expand Up @@ -972,6 +973,8 @@ list (APPEND PUBLIC_HEADER_FILES
opm/simulators/wells/BlackoilWellModelGasLift_impl.hpp
opm/simulators/wells/BlackoilWellModelGeneric.hpp
opm/simulators/wells/BlackoilWellModelGuideRates.hpp
opm/simulators/wells/BlackoilWellModelNldd.hpp
opm/simulators/wells/BlackoilWellModelNldd_impl.hpp
opm/simulators/wells/BlackoilWellModelRestart.hpp
opm/simulators/wells/BlackoilWellModelWBP.hpp
opm/simulators/wells/ConnectionIndexMap.hpp
Expand Down
39 changes: 24 additions & 15 deletions opm/simulators/flow/BlackoilModelNldd.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -50,6 +50,8 @@

#include <opm/simulators/utils/ComponentName.hpp>

#include <opm/simulators/wells/BlackoilWellModelNldd.hpp>

#include <fmt/format.h>

#include <algorithm>
Expand All @@ -64,6 +66,7 @@
#include <ios>
#include <memory>
#include <numeric>
#include <set>
#include <sstream>
#include <string>
#include <type_traits>
Expand All @@ -76,7 +79,8 @@ template<class TypeTag> class BlackoilModel;

/// A NLDD implementation for three-phase black oil.
template <class TypeTag>
class BlackoilModelNldd {
class BlackoilModelNldd
{
public:
using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
Expand All @@ -98,11 +102,16 @@ class BlackoilModelNldd {
//! \param param param Model parameters
//! \param compNames Names of the solution components
explicit BlackoilModelNldd(BlackoilModel<TypeTag>& model)
: model_(model), rank_(model_.simulator().vanguard().grid().comm().rank())
: model_(model)
, wellModel_(model.wellModel())
, rank_(model_.simulator().vanguard().grid().comm().rank())
{
// Create partitions.
const auto& [partition_vector, num_domains] = this->partitionCells();

// Set nldd handler in main well model
model.wellModel().setNlddAdapter(&wellModel_);

// Scan through partitioning to get correct size for each.
std::vector<int> sizes(num_domains, 0);
for (const auto& p : partition_vector) {
Expand Down Expand Up @@ -179,7 +188,7 @@ class BlackoilModelNldd {
void prepareStep()
{
// Setup domain->well mapping.
model_.wellModel().setupDomains(domains_);
wellModel_.setupDomains(domains_);
}

//! \brief Do one non-linear NLDD iteration.
Expand Down Expand Up @@ -352,7 +361,6 @@ class BlackoilModelNldd {
}

private:

//! \brief Solve the equation system for a single domain.
std::pair<SimulatorReportSingle, ConvergenceReport>
solveDomain(const Domain& domain,
Expand Down Expand Up @@ -380,9 +388,9 @@ class BlackoilModelNldd {
modelSimulator.model().newtonMethod().setIterationIndex(iter);
// TODO: we should have a beginIterationLocal function()
// only handling the well model for now
modelSimulator.problem().wellModel().assembleDomain(modelSimulator.model().newtonMethod().numIterations(),
modelSimulator.timeStepSize(),
domain);
wellModel_.assemble(modelSimulator.model().newtonMethod().numIterations(),
modelSimulator.timeStepSize(),
domain);
// Assemble reservoir locally.
report += this->assembleReservoirDomain(domain);
report.assemble_time += detailTimer.stop();
Expand Down Expand Up @@ -445,9 +453,9 @@ class BlackoilModelNldd {
// TODO: we should have a beginIterationLocal function()
// only handling the well model for now
// Assemble reservoir locally.
modelSimulator.problem().wellModel().assembleDomain(modelSimulator.model().newtonMethod().numIterations(),
modelSimulator.timeStepSize(),
domain);
wellModel_.assemble(modelSimulator.model().newtonMethod().numIterations(),
modelSimulator.timeStepSize(),
domain);
report += this->assembleReservoirDomain(domain);
report.assemble_time += detailTimer.stop();

Expand Down Expand Up @@ -736,7 +744,7 @@ class BlackoilModelNldd {
logger,
B_avg,
residual_norms);
report += model_.wellModel().getDomainWellConvergence(domain, B_avg, logger);
report += wellModel_.getWellConvergence(domain, B_avg, logger);
return report;
}

Expand Down Expand Up @@ -815,7 +823,7 @@ class BlackoilModelNldd {
const SimulatorTimerInterface& timer,
const Domain& domain)
{
auto initial_local_well_primary_vars = model_.wellModel().getPrimaryVarsDomain(domain.index);
auto initial_local_well_primary_vars = wellModel_.getPrimaryVarsDomain(domain.index);
auto initial_local_solution = Details::extractVector(solution, domain.cells);
auto res = solveDomain(domain, timer, logger, iteration, false);
local_report = res.first;
Expand All @@ -825,7 +833,7 @@ class BlackoilModelNldd {
Details::setGlobal(initial_local_solution, domain.cells, solution);
model_.simulator().model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0, domain);
} else {
model_.wellModel().setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
wellModel_.setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
Details::setGlobal(initial_local_solution, domain.cells, solution);
model_.simulator().model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0, domain);
}
Expand All @@ -840,7 +848,7 @@ class BlackoilModelNldd {
const SimulatorTimerInterface& timer,
const Domain& domain)
{
auto initial_local_well_primary_vars = model_.wellModel().getPrimaryVarsDomain(domain.index);
auto initial_local_well_primary_vars = wellModel_.getPrimaryVarsDomain(domain.index);
auto initial_local_solution = Details::extractVector(solution, domain.cells);
auto res = solveDomain(domain, timer, logger, iteration, true);
local_report = res.first;
Expand Down Expand Up @@ -881,7 +889,7 @@ class BlackoilModelNldd {
auto local_solution = Details::extractVector(solution, domain.cells);
Details::setGlobal(local_solution, domain.cells, locally_solved);
} else {
model_.wellModel().setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
wellModel_.setPrimaryVarsDomain(domain.index, initial_local_well_primary_vars);
Details::setGlobal(initial_local_solution, domain.cells, solution);
model_.simulator().model().invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0, domain);
}
Expand Down Expand Up @@ -996,6 +1004,7 @@ class BlackoilModelNldd {
}

BlackoilModel<TypeTag>& model_; //!< Reference to model
BlackoilWellModelNldd<TypeTag> wellModel_; //!< NLDD well model adapter
std::vector<Domain> domains_; //!< Vector of subdomains
std::vector<std::unique_ptr<Mat>> domain_matrices_; //!< Vector of matrix operator for each subdomain
std::vector<ISTLSolverType> domain_linsolvers_; //!< Vector of linear solvers for each domain
Expand Down
16 changes: 13 additions & 3 deletions opm/simulators/flow/SubDomain.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -57,8 +57,7 @@ namespace Opm

/// Representing a part of a grid, in a way suitable for performing
/// local solves.
template <class Grid>
struct SubDomain
struct SubDomainIndices
{
// The index of a subdomain is arbitrary, but can be used by the
// solvers to keep track of well locations etc.
Expand All @@ -72,10 +71,21 @@ namespace Opm
std::vector<bool> interior;
// Enables subdomain solves and linearization using the generic linearization
// approach (i.e. FvBaseLinearizer as opposed to TpfaLinearizer).
SubDomainIndices(const int i, std::vector<int>&& c, std::vector<bool>&& in)
: index(i), cells(std::move(c)), interior(std::move(in))
{}
};

/// Representing a part of a grid, in a way suitable for performing
/// local solves.
template <class Grid>
struct SubDomain : public SubDomainIndices
{
Dune::SubGridPart<Grid> view;
// Constructor that moves from its argument.
SubDomain(const int i, std::vector<int>&& c, std::vector<bool>&& in, Dune::SubGridPart<Grid>&& v)
: index(i), cells(std::move(c)), interior(std::move(in)), view(std::move(v))
: SubDomainIndices(i, std::move(c), std::move(in))
, view(std::move(v))
{}
};

Expand Down
8 changes: 6 additions & 2 deletions opm/simulators/linalg/WellOperators.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -188,7 +188,8 @@ class DomainWellModelAsLinearOperator : public WellModelAsLinearOperator<WellMod
std::size_t well_index = 0;
for (const auto& well : this->wellMod_) {
if (this->wellMod_.well_domain().at(well->name()) == domainIndex_) {
this->applySingleWell(x, y, well, this->wellMod_.well_local_cells()[well_index]);
this->applySingleWell(x, y, well,
this->wellMod_.well_local_cells()[well_index]);
}
++well_index;
}
Expand All @@ -199,7 +200,10 @@ class DomainWellModelAsLinearOperator : public WellModelAsLinearOperator<WellMod
const bool use_well_weights) const override
{
OPM_TIMEBLOCK(addWellPressureEquations);
this->wellMod_.addWellPressureEquationsDomain(jacobian, weights, use_well_weights, domainIndex_);
this->wellMod_.addWellPressureEquationsDomain(jacobian,
weights,
use_well_weights,
domainIndex_);
}

private:
Expand Down
87 changes: 49 additions & 38 deletions opm/simulators/wells/BlackoilWellModel.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -30,8 +30,6 @@

#include <opm/common/OpmLog/OpmLog.hpp>

#include <opm/grid/utility/SparseTable.hpp>

#include <opm/input/eclipse/Schedule/Group/Group.hpp>
#include <opm/input/eclipse/Schedule/Group/GuideRate.hpp>
#include <opm/input/eclipse/Schedule/Schedule.hpp>
Expand All @@ -41,7 +39,6 @@

#include <opm/simulators/flow/countGlobalCells.hpp>
#include <opm/simulators/flow/FlowBaseVanguard.hpp>
#include <opm/simulators/flow/SubDomain.hpp>

#include <opm/simulators/linalg/matrixblock.hh>

Expand Down Expand Up @@ -81,6 +78,9 @@

namespace Opm {

template<class Scalar> class BlackoilWellModelNldd;
template<class T> class SparseTable;

#if COMPILE_GPU_BRIDGE
template<class Scalar> class WellContributions;
#endif
Expand Down Expand Up @@ -132,8 +132,6 @@ template<class Scalar> class WellContributions;
using AverageRegionalPressureType = RegionAverageCalculator::
AverageRegionalPressure<FluidSystem, std::vector<int> >;

using Domain = SubDomain<Grid>;

explicit BlackoilWellModel(Simulator& simulator);

void init();
Expand Down Expand Up @@ -248,11 +246,6 @@ template<class Scalar> class WellContributions;
// Check if well equations is converged.
ConvergenceReport getWellConvergence(const std::vector<Scalar>& B_avg, const bool checkWellGroupControls = false) const;

// Check if well equations are converged locally.
ConvergenceReport getDomainWellConvergence(const Domain& domain,
const std::vector<Scalar>& B_avg,
DeferredLogger& local_deferredLogger) const;

const SimulatorReportSingle& lastReport() const;

void addWellContributions(SparseMatrixAdapter& jacobian) const;
Expand All @@ -274,7 +267,6 @@ template<class Scalar> class WellContributions;
// at the beginning of each time step (Not report step)
void prepareTimeStep(DeferredLogger& deferred_logger);
void initPrimaryVariablesEvaluation() const;
void initPrimaryVariablesEvaluationDomain(const Domain& domain) const;

std::pair<bool, bool>
updateWellControls(const bool mandatory_network_balance, DeferredLogger& deferred_logger, const bool relax_network_tolerance = false);
Expand All @@ -292,40 +284,56 @@ template<class Scalar> class WellContributions;

using PressureMatrix = Dune::BCRSMatrix<Opm::MatrixBlock<Scalar, 1, 1>>;

void addWellPressureEquations(PressureMatrix& jacobian, const BVector& weights,const bool use_well_weights) const;

void addWellPressureEquationsDomain([[maybe_unused]] PressureMatrix& jacobian,
[[maybe_unused]] const BVector& weights,
[[maybe_unused]] const bool use_well_weights,
[[maybe_unused]] const int domainIndex) const;


void addWellPressureEquations(PressureMatrix& jacobian,
const BVector& weights,
const bool use_well_weights) const;
void addWellPressureEquationsStruct(PressureMatrix& jacobian) const;
void addWellPressureEquationsDomain(PressureMatrix& jacobian,
const BVector& weights,
const bool use_well_weights,
const int domainIndex) const
{
if (!nldd_) {
OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
}
return nldd_->addWellPressureEquations(jacobian,
weights,
use_well_weights,
domainIndex);
}

/// \brief Get list of local nonshut wells
const std::vector<WellInterfacePtr>& localNonshutWells() const
{
return well_container_;
}

// prototype for assemble function for ASPIN solveLocal()
// will try to merge back to assemble() when done prototyping
void assembleDomain(const int iterationIdx,
const double dt,
const Domain& domain);
void updateWellControlsDomain(DeferredLogger& deferred_logger, const Domain& domain);

void setupDomains(const std::vector<Domain>& domains);

const SparseTable<int>& well_local_cells() const
{
return well_local_cells_;
if (!nldd_) {
OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
}
return nldd_->well_local_cells();
}

const std::map<std::string, int>& well_domain() const
{
if (!nldd_) {
OPM_THROW(std::logic_error, "Attempt to access NLDD data without a NLDD solver");
}

return nldd_->well_domain();
}

auto begin() const { return well_container_.begin(); }
auto end() const { return well_container_.end(); }
bool empty() const { return well_container_.empty(); }

bool addMatrixContributions() const { return param_.matrix_add_well_contributions_; }
bool addMatrixContributions() const
{ return param_.matrix_add_well_contributions_; }

int numStrictIterations() const
{ return param_.strict_outer_iter_wells_; }

int compressedIndexForInterior(int cartesian_cell_idx) const override
{
Expand All @@ -339,7 +347,16 @@ template<class Scalar> class WellContributions;
// using the solution x to recover the solution xw for wells and applying
// xw to update Well State
void recoverWellSolutionAndUpdateWellStateDomain(const BVector& x,
const Domain& domain);
const int domainIdx);

const Grid& grid() const
{ return simulator_.vanguard().grid(); }

const Simulator& simulator() const
{ return simulator_; }

void setNlddAdapter(BlackoilWellModelNldd<TypeTag>* mod)
{ nldd_ = mod; }

protected:
Simulator& simulator_;
Expand Down Expand Up @@ -385,12 +402,6 @@ template<class Scalar> class WellContributions;

std::vector<Scalar> B_avg_{};

// Store the local index of the wells perforated cells in the domain, if using subdomains
SparseTable<int> well_local_cells_;

const Grid& grid() const
{ return simulator_.vanguard().grid(); }

const EquilGrid& equilGrid() const
{ return simulator_.vanguard().equilGrid(); }

Expand Down Expand Up @@ -464,7 +475,6 @@ template<class Scalar> class WellContributions;
int reportStepIndex() const;

void assembleWellEq(const double dt, DeferredLogger& deferred_logger);
void assembleWellEqDomain(const double dt, const Domain& domain, DeferredLogger& deferred_logger);

void prepareWellsBeforeAssembling(const double dt, DeferredLogger& deferred_logger);

Expand Down Expand Up @@ -495,6 +505,7 @@ template<class Scalar> class WellContributions;
BlackoilWellModel(Simulator& simulator, const PhaseUsage& pu);

BlackoilWellModelGasLift<TypeTag> gaslift_;
BlackoilWellModelNldd<TypeTag>* nldd_ = nullptr; //!< NLDD well model adapter (not owned)

// These members are used to avoid reallocation in specific functions
// instead of using local variables.
Expand Down
Loading