2016-04-06 08:23:43 -05:00
|
|
|
/*
|
|
|
|
Copyright 2016 SINTEF ICT, Applied Mathematics.
|
2016-04-06 08:31:58 -05:00
|
|
|
Copyright 2016 Statoil ASA.
|
2016-04-06 08:23:43 -05:00
|
|
|
|
|
|
|
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_STANDARDWELLS_HEADER_INCLUDED
|
|
|
|
#define OPM_STANDARDWELLS_HEADER_INCLUDED
|
|
|
|
|
2016-09-23 04:50:11 -05:00
|
|
|
#include <dune/common/parallel/mpihelper.hh>
|
|
|
|
|
2016-06-29 04:43:11 -05:00
|
|
|
#include <opm/common/OpmLog/OpmLog.hpp>
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
#include <opm/common/utility/platform_dependent/disable_warnings.h>
|
|
|
|
#include <Eigen/Eigen>
|
|
|
|
#include <Eigen/Sparse>
|
|
|
|
#include <opm/common/utility/platform_dependent/reenable_warnings.h>
|
|
|
|
|
|
|
|
#include <cassert>
|
2016-07-05 06:02:06 -05:00
|
|
|
#include <tuple>
|
2016-04-06 08:23:43 -05:00
|
|
|
|
2016-06-27 09:14:28 -05:00
|
|
|
#include <opm/parser/eclipse/EclipseState/Schedule/Schedule.hpp>
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
#include <opm/core/wells.h>
|
2016-06-27 09:14:28 -05:00
|
|
|
#include <opm/core/wells/DynamicListEconLimited.hpp>
|
2016-09-30 05:57:24 -05:00
|
|
|
#include <opm/core/wells/WellCollection.hpp>
|
2016-04-06 08:23:43 -05:00
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
2016-04-07 09:17:25 -05:00
|
|
|
#include <opm/autodiff/AutoDiffHelpers.hpp>
|
2016-12-29 09:35:24 -06:00
|
|
|
#include <opm/autodiff/BlackoilPropsAdFromDeck.hpp>
|
2016-09-23 04:50:11 -05:00
|
|
|
#include <opm/simulators/WellSwitchingLogger.hpp>
|
2016-04-06 08:23:43 -05:00
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
2016-04-07 09:17:25 -05:00
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
/// Class for handling the standard well model.
|
|
|
|
class StandardWells {
|
2016-06-20 04:02:49 -05:00
|
|
|
public:
|
2016-04-06 08:23:43 -05:00
|
|
|
struct WellOps {
|
|
|
|
explicit WellOps(const Wells* wells);
|
|
|
|
Eigen::SparseMatrix<double> w2p; // well -> perf (scatter)
|
|
|
|
Eigen::SparseMatrix<double> p2w; // perf -> well (gather)
|
|
|
|
std::vector<int> well_cells; // the set of perforated cells
|
|
|
|
};
|
|
|
|
|
2016-04-12 04:35:03 -05:00
|
|
|
// --------- Types ---------
|
|
|
|
using ADB = AutoDiffBlock<double>;
|
|
|
|
using Vector = ADB::V;
|
2016-09-23 04:50:11 -05:00
|
|
|
using Communication =
|
|
|
|
Dune::CollectiveCommunication<typename Dune::MPIHelper
|
|
|
|
::MPICommunicator>;
|
2016-04-12 04:35:03 -05:00
|
|
|
|
|
|
|
// copied from BlackoilModelBase
|
|
|
|
// should put to somewhere better
|
|
|
|
using DataBlock = Eigen::Array<double,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::RowMajor>;
|
2016-04-06 08:23:43 -05:00
|
|
|
// --------- Public methods ---------
|
2016-09-30 10:53:55 -05:00
|
|
|
StandardWells(const Wells* wells_arg, WellCollection* well_collection);
|
2016-05-02 09:21:14 -05:00
|
|
|
|
2016-12-29 09:35:24 -06:00
|
|
|
void init(const BlackoilPropsAdFromDeck* fluid_arg,
|
2016-05-02 09:21:14 -05:00
|
|
|
const std::vector<bool>* active_arg,
|
2016-05-09 07:00:21 -05:00
|
|
|
const std::vector<PhasePresence>* pc_arg,
|
|
|
|
const VFPProperties* vfp_properties_arg,
|
|
|
|
const double gravity_arg,
|
2016-05-09 08:27:09 -05:00
|
|
|
const Vector& depth_arg);
|
2016-04-28 09:20:18 -05:00
|
|
|
|
|
|
|
const WellOps& wellOps() const;
|
|
|
|
|
2016-06-27 05:14:17 -05:00
|
|
|
int numPhases() const { return wells().number_of_phases; };
|
2016-04-06 08:23:43 -05:00
|
|
|
|
|
|
|
const Wells& wells() const;
|
|
|
|
|
2016-06-27 06:03:30 -05:00
|
|
|
const Wells* wellsPointer() const;
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
/// return true if wells are available in the reservoir
|
|
|
|
bool wellsActive() const;
|
|
|
|
void setWellsActive(const bool wells_active);
|
|
|
|
/// return true if wells are available on this process
|
|
|
|
bool localWellsActive() const;
|
|
|
|
|
2016-05-11 03:35:11 -05:00
|
|
|
int numWellVars() const;
|
2016-04-06 08:23:43 -05:00
|
|
|
|
|
|
|
/// Density of each well perforation
|
2016-04-12 04:35:03 -05:00
|
|
|
Vector& wellPerforationDensities(); // mutable version kept for BlackoilMultisegmentModel
|
2016-04-06 08:23:43 -05:00
|
|
|
const Vector& wellPerforationDensities() const;
|
|
|
|
|
|
|
|
/// Diff to bhp for each well perforation.
|
2016-04-12 04:35:03 -05:00
|
|
|
Vector& wellPerforationPressureDiffs(); // mutable version kept for BlackoilMultisegmentModel
|
2016-04-06 08:23:43 -05:00
|
|
|
const Vector& wellPerforationPressureDiffs() const;
|
|
|
|
|
2016-04-11 10:10:15 -05:00
|
|
|
template <class ReservoirResidualQuant, class SolutionState>
|
|
|
|
void extractWellPerfProperties(const SolutionState& state,
|
|
|
|
const std::vector<ReservoirResidualQuant>& rq,
|
2016-04-08 09:48:31 -05:00
|
|
|
std::vector<ADB>& mob_perfcells,
|
|
|
|
std::vector<ADB>& b_perfcells) const;
|
|
|
|
|
2016-04-08 10:26:07 -05:00
|
|
|
template <class SolutionState>
|
|
|
|
void computeWellFlux(const SolutionState& state,
|
|
|
|
const std::vector<ADB>& mob_perfcells,
|
|
|
|
const std::vector<ADB>& b_perfcells,
|
|
|
|
Vector& aliveWells,
|
|
|
|
std::vector<ADB>& cq_s) const;
|
|
|
|
|
2016-04-11 03:46:21 -05:00
|
|
|
template <class SolutionState, class WellState>
|
|
|
|
void updatePerfPhaseRatesAndPressures(const std::vector<ADB>& cq_s,
|
|
|
|
const SolutionState& state,
|
|
|
|
WellState& xw) const;
|
|
|
|
|
2016-04-11 08:45:03 -05:00
|
|
|
template <class WellState>
|
|
|
|
void updateWellState(const Vector& dwells,
|
|
|
|
const double dpmaxrel,
|
|
|
|
WellState& well_state);
|
|
|
|
|
2016-04-19 10:09:50 -05:00
|
|
|
template <class WellState>
|
2016-06-07 11:56:37 -05:00
|
|
|
void updateWellControls(WellState& xw) const;
|
2016-04-19 10:09:50 -05:00
|
|
|
|
|
|
|
// TODO: should LinearisedBlackoilResidual also be a template class?
|
|
|
|
template <class SolutionState>
|
|
|
|
void addWellFluxEq(const std::vector<ADB>& cq_s,
|
|
|
|
const SolutionState& state,
|
|
|
|
LinearisedBlackoilResidual& residual);
|
|
|
|
|
|
|
|
// TODO: some parameters, like gravity, maybe it is better to put in the member list
|
|
|
|
template <class SolutionState, class WellState>
|
|
|
|
void addWellControlEq(const SolutionState& state,
|
|
|
|
const WellState& xw,
|
|
|
|
const Vector& aliveWells,
|
|
|
|
LinearisedBlackoilResidual& residual);
|
|
|
|
|
|
|
|
template <class SolutionState, class WellState>
|
|
|
|
void computeWellConnectionPressures(const SolutionState& state,
|
2016-05-09 08:27:09 -05:00
|
|
|
const WellState& xw);
|
2016-04-19 10:09:50 -05:00
|
|
|
|
2016-04-21 08:46:20 -05:00
|
|
|
// state0 is non-constant, while it will not be used outside of the function
|
|
|
|
template <class SolutionState, class WellState>
|
|
|
|
void
|
2016-05-11 06:13:08 -05:00
|
|
|
computeWellPotentials(const std::vector<ADB>& mob_perfcells,
|
2016-04-21 08:46:20 -05:00
|
|
|
const std::vector<ADB>& b_perfcells,
|
2017-04-05 05:52:07 -05:00
|
|
|
const WellState& well_state,
|
2016-05-11 06:13:08 -05:00
|
|
|
SolutionState& state0,
|
2017-04-05 05:52:07 -05:00
|
|
|
std::vector<double>& well_potentials) const;
|
2016-04-19 04:20:18 -05:00
|
|
|
|
2016-05-12 08:14:20 -05:00
|
|
|
template <class SolutionState>
|
|
|
|
void
|
|
|
|
variableStateExtractWellsVars(const std::vector<int>& indices,
|
|
|
|
std::vector<ADB>& vars,
|
|
|
|
SolutionState& state) const;
|
2016-04-18 04:22:19 -05:00
|
|
|
|
2016-04-21 10:42:50 -05:00
|
|
|
void
|
|
|
|
variableStateWellIndices(std::vector<int>& indices,
|
|
|
|
int& next) const;
|
|
|
|
|
2016-04-21 11:03:16 -05:00
|
|
|
std::vector<int>
|
|
|
|
variableWellStateIndices() const;
|
2016-04-07 09:17:25 -05:00
|
|
|
|
2016-04-22 03:51:17 -05:00
|
|
|
template <class WellState>
|
|
|
|
void
|
|
|
|
variableWellStateInitials(const WellState& xw,
|
|
|
|
std::vector<Vector>& vars0) const;
|
|
|
|
|
2016-06-20 04:02:49 -05:00
|
|
|
/// If set, computeWellFlux() will additionally store the
|
|
|
|
/// total reservoir volume perforation fluxes.
|
|
|
|
void setStoreWellPerforationFluxesFlag(const bool store_fluxes);
|
|
|
|
|
|
|
|
/// Retrieves the stored fluxes. It is an error to call this
|
|
|
|
/// unless setStoreWellPerforationFluxesFlag(true) has been
|
|
|
|
/// called.
|
|
|
|
const Vector& getStoredWellPerforationFluxes() const;
|
|
|
|
|
2016-06-27 09:14:28 -05:00
|
|
|
/// upate the dynamic lists related to economic limits
|
|
|
|
template<class WellState>
|
|
|
|
void
|
2016-10-14 02:23:26 -05:00
|
|
|
updateListEconLimited(const Schedule& schedule,
|
2016-06-27 09:14:28 -05:00
|
|
|
const int current_step,
|
|
|
|
const Wells* wells,
|
|
|
|
const WellState& well_state,
|
|
|
|
DynamicListEconLimited& list_econ_limited) const;
|
|
|
|
|
2016-10-06 11:00:36 -05:00
|
|
|
|
2016-10-13 11:16:19 -05:00
|
|
|
WellCollection* wellCollection() const;
|
2016-10-06 11:00:36 -05:00
|
|
|
|
2016-10-19 03:40:11 -05:00
|
|
|
void calculateEfficiencyFactors();
|
|
|
|
|
|
|
|
const Vector& wellPerfEfficiencyFactors() const;
|
2016-10-18 10:32:42 -05:00
|
|
|
|
2017-11-08 06:57:36 -06:00
|
|
|
// just return the passed well state
|
|
|
|
template<class WellState>
|
|
|
|
const WellState& wellState(const WellState& well_state) const { return well_state; }
|
|
|
|
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
protected:
|
|
|
|
bool wells_active_;
|
|
|
|
const Wells* wells_;
|
|
|
|
const WellOps wops_;
|
2016-10-13 11:16:19 -05:00
|
|
|
// It will probably need to be updated during running time.
|
2016-09-30 10:53:55 -05:00
|
|
|
WellCollection* well_collection_;
|
2016-05-02 09:21:14 -05:00
|
|
|
|
2016-10-18 10:32:42 -05:00
|
|
|
// The efficiency factor for each connection
|
|
|
|
// It is specified based on wells and groups
|
|
|
|
// By default, they should all be one.
|
2016-10-19 03:40:11 -05:00
|
|
|
Vector well_perforation_efficiency_factors_;
|
2016-10-18 10:32:42 -05:00
|
|
|
|
2016-12-29 09:35:24 -06:00
|
|
|
const BlackoilPropsAdFromDeck* fluid_;
|
2016-05-02 09:21:14 -05:00
|
|
|
const std::vector<bool>* active_;
|
|
|
|
const std::vector<PhasePresence>* phase_condition_;
|
2016-05-09 07:00:21 -05:00
|
|
|
const VFPProperties* vfp_properties_;
|
|
|
|
double gravity_;
|
2016-05-09 08:27:09 -05:00
|
|
|
// the depth of the all the cell centers
|
|
|
|
// for standard Wells, it the same with the perforation depth
|
|
|
|
Vector perf_cell_depth_;
|
2016-05-02 09:21:14 -05:00
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
Vector well_perforation_densities_;
|
|
|
|
Vector well_perforation_pressure_diffs_;
|
2016-04-29 05:51:32 -05:00
|
|
|
|
2016-06-20 04:02:49 -05:00
|
|
|
bool store_well_perforation_fluxes_;
|
|
|
|
Vector well_perforation_fluxes_;
|
|
|
|
|
2016-04-29 05:51:32 -05:00
|
|
|
// protected methods
|
|
|
|
template <class SolutionState, class WellState>
|
|
|
|
void computePropertiesForWellConnectionPressures(const SolutionState& state,
|
|
|
|
const WellState& xw,
|
|
|
|
std::vector<double>& b_perf,
|
|
|
|
std::vector<double>& rsmax_perf,
|
|
|
|
std::vector<double>& rvmax_perf,
|
|
|
|
std::vector<double>& surf_dens_perf);
|
|
|
|
|
|
|
|
template <class WellState>
|
|
|
|
void computeWellConnectionDensitesPressures(const WellState& xw,
|
|
|
|
const std::vector<double>& b_perf,
|
|
|
|
const std::vector<double>& rsmax_perf,
|
|
|
|
const std::vector<double>& rvmax_perf,
|
|
|
|
const std::vector<double>& surf_dens_perf,
|
|
|
|
const std::vector<double>& depth_perf,
|
|
|
|
const double grav);
|
|
|
|
|
2016-06-28 09:13:08 -05:00
|
|
|
|
|
|
|
template <class WellState>
|
|
|
|
bool checkRateEconLimits(const WellEconProductionLimits& econ_production_limits,
|
|
|
|
const WellState& well_state,
|
|
|
|
const int well_number) const;
|
|
|
|
|
2016-06-28 10:27:25 -05:00
|
|
|
using WellMapType = typename WellState::WellMapType;
|
2016-07-05 08:28:02 -05:00
|
|
|
using WellMapEntryType = typename WellState::mapentry_t;
|
2016-06-28 10:27:25 -05:00
|
|
|
|
2016-07-05 06:02:06 -05:00
|
|
|
// a tuple type for ratio limit check.
|
|
|
|
// first value indicates whether ratio limit is violated, when the ratio limit is not violated, the following three
|
|
|
|
// values should not be used.
|
|
|
|
// second value indicates whehter there is only one connection left.
|
|
|
|
// third value indicates the indx of the worst-offending connection.
|
|
|
|
// the last value indicates the extent of the violation for the worst-offending connection, which is defined by
|
|
|
|
// the ratio of the actual value to the value of the violated limit.
|
|
|
|
using RatioCheckTuple = std::tuple<bool, bool, int, double>;
|
|
|
|
|
|
|
|
enum ConnectionIndex {
|
|
|
|
INVALIDCONNECTION = -10000
|
|
|
|
};
|
|
|
|
|
2016-06-28 10:27:25 -05:00
|
|
|
|
2016-06-28 10:59:50 -05:00
|
|
|
template <class WellState>
|
2016-07-05 06:02:06 -05:00
|
|
|
RatioCheckTuple checkRatioEconLimits(const WellEconProductionLimits& econ_production_limits,
|
|
|
|
const WellState& well_state,
|
2016-07-05 08:28:02 -05:00
|
|
|
const WellMapEntryType& map_entry) const;
|
2016-06-28 10:59:50 -05:00
|
|
|
|
2016-06-28 10:27:25 -05:00
|
|
|
template <class WellState>
|
2016-07-05 06:02:06 -05:00
|
|
|
RatioCheckTuple checkMaxWaterCutLimit(const WellEconProductionLimits& econ_production_limits,
|
|
|
|
const WellState& well_state,
|
2016-07-05 08:28:02 -05:00
|
|
|
const WellMapEntryType& map_entry) const;
|
2016-06-28 10:27:25 -05:00
|
|
|
|
2016-11-09 08:27:20 -06:00
|
|
|
template <class WellState>
|
|
|
|
void updateWellStateWithTarget(const WellControls* wc,
|
|
|
|
const int current,
|
|
|
|
const int well_index,
|
|
|
|
WellState& xw) const;
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
} // namespace Opm
|
2016-04-07 08:41:08 -05:00
|
|
|
|
|
|
|
#include "StandardWells_impl.hpp"
|
|
|
|
|
2016-04-06 08:23:43 -05:00
|
|
|
#endif
|