Files
opm-simulators/opm/autodiff/WellInterface.hpp
Kai Bao ab67635134 adding applySolutionWellState to apply solution from reservoir
to update well state.

With this way, the BlackoilModelEbos does not need to know the data type
assocated with different well type.

It is not well tested yet.
2017-08-25 14:09:26 +02:00

265 lines
8.9 KiB
C++

/*
Copyright 2017 SINTEF ICT, Applied Mathematics.
Copyright 2017 Statoil 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_WELLINTERFACE_HEADER_INCLUDED
#define OPM_WELLINTERFACE_HEADER_INCLUDED
#include <opm/common/OpmLog/OpmLog.hpp>
#include <opm/parser/eclipse/EclipseState/Schedule/Well.hpp>
#include <opm/core/wells.h>
#include <opm/core/well_controls.h>
#include <opm/core/props/BlackoilPhases.hpp>
#include <opm/autodiff/VFPProperties.hpp>
#include <opm/autodiff/VFPInjProperties.hpp>
#include <opm/autodiff/VFPProdProperties.hpp>
#include <opm/autodiff/WellHelpers.hpp>
#include <opm/autodiff/WellStateFullyImplicitBlackoilDense.hpp>
#include <opm/autodiff/BlackoilModelParameters.hpp>
#include <opm/simulators/WellSwitchingLogger.hpp>
#include<dune/common/fmatrix.hh>
#include<dune/istl/bcrsmatrix.hh>
#include<dune/istl/matrixmatrix.hh>
#include <opm/material/densead/Math.hpp>
#include <opm/material/densead/Evaluation.hpp>
#include <string>
#include <memory>
#include <vector>
#include <cassert>
namespace Opm
{
template<typename TypeTag>
class WellInterface
{
public:
using WellState = WellStateFullyImplicitBlackoilDense;
typedef BlackoilModelParameters ModelParameters;
typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
typedef typename GET_PROP_TYPE(TypeTag, Indices) BlackoilIndices;
typedef typename GET_PROP_TYPE(TypeTag, IntensiveQuantities) IntensiveQuantities;
typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
static const int numEq = BlackoilIndices::numEq;
typedef double Scalar;
typedef Dune::FieldVector<Scalar, numEq > VectorBlockType;
typedef Dune::FieldMatrix<Scalar, numEq, numEq > MatrixBlockType;
typedef Dune::BCRSMatrix <MatrixBlockType> Mat;
typedef Dune::BlockVector<VectorBlockType> BVector;
typedef DenseAd::Evaluation<double, /*size=*/numEq> Eval;
static const bool has_solvent = GET_PROP_VALUE(TypeTag, EnableSolvent);
/// Constructor
// TODO: Wel can be reference.
WellInterface(const Well* well, const int time_step, const Wells* wells);
/// Well name.
const std::string& name() const;
/// The index of the well in Wells struct
// It is used to locate the inforation in Wells and also WellState for now.
int indexOfWell() const;
/// Well type, INJECTOR or PRODUCER.
WellType wellType() const;
/// number of phases
int numberOfPhases() const;
/// Component fractions for each phase for the well
const std::vector<double>& compFrac() const;
/// Well controls
// TODO: later to see whether we need to return const.
WellControls* wellControls() const;
/// Number of the perforations
int numberOfPerforations() const;
/// Well productivity index for each perforation.
const std::vector<double>& wellIndex() const;
/// Depth of perforations
const std::vector<double>& perfDepth() const;
/// Indices of the grid cells/blocks that perforations are completed within.
const std::vector<int>& wellCells() const;
// TODO: the following function should be able to be removed by refactoring the well class
// It is probably only needed for StandardWell
/// the densities of the fluid in each perforation
virtual const std::vector<double>& perfDensities() const = 0;
virtual std::vector<double>& perfDensities() = 0;
/// the pressure difference between different perforations
virtual const std::vector<double>& perfPressureDiffs() const = 0;
virtual std::vector<double>& perfPressureDiffs() = 0;
// TODO: the parameters need to be optimized/adjusted
virtual void init(const PhaseUsage* phase_usage_arg,
const std::vector<bool>* active_arg,
const VFPProperties* vfp_properties_arg,
const double gravity_arg,
const int num_cells);
// TODO: temporary
virtual void setWellVariables(const WellState& well_state) = 0;
const std::vector<bool>& active() const;
const PhaseUsage& phaseUsage() const;
int flowPhaseToEbosCompIdx( const int phaseIdx ) const;
int flowToEbosPvIdx( const int flowPv ) const;
int flowPhaseToEbosPhaseIdx( const int phaseIdx ) const;
int numPhases() const;
int numComponents() const;
// simply returning allow_cf_
// TODO: to check whether needed, it causes name problem with the crossFlowAllowed
bool allowCrossFlow() const;
virtual bool crossFlowAllowed(const Simulator& ebosSimulator) const = 0;
// TODO: for this kind of function, maybe can make a function with parameter perf
const std::vector<int>& saturationTableNumber() const;
const double wsolvent() const;
virtual bool getWellConvergence(Simulator& ebosSimulator,
const std::vector<double>& B_avg,
const ModelParameters& param) const = 0;
virtual void wellEqIteration(Simulator& ebosSimulator,
const ModelParameters& param,
WellState& well_state) = 0;
virtual void assembleWellEq(Simulator& ebosSimulator,
const double dt,
WellState& well_state,
bool only_wells) = 0;
virtual void updateWellStateWithTarget(const int current,
WellState& xw) const = 0;
virtual void updateWellControl(WellState& xw) const = 0;
virtual void computeAccumWell() = 0;
// TODO: it should come with a different name
// for MS well, the definition is different and should not use this name anymore
virtual void computeWellConnectionPressures(const Simulator& ebosSimulator,
const WellState& xw) = 0;
// Ax = Ax - C D^-1 B x
virtual void apply(const BVector& x, BVector& Ax) const = 0;
// r = r - C D^-1 Rw
virtual void apply(BVector& r) const = 0;
// using the solution x to recover the solution xw for wells and applying
// xw to update Well State
virtual void applySolutionWellState(const BVector& x, const ModelParameters& param,
WellState& well_state) const = 0;
protected:
// TODO: some variables shared by all the wells should be made static
// well name
std::string name_;
// the index of well in Wells struct
int index_of_well_;
// well type
// INJECTOR or PRODUCER
enum WellType well_type_;
// whether the well allows crossflow
bool allow_cf_;
// number of phases
int number_of_phases_;
// component fractions for each well
// typically, it should apply to injection wells
std::vector<double> comp_frac_;
// controls for this well
// TODO: later will check whehter to let it stay with pointer
struct WellControls* well_controls_;
// number of the perforations for this well
int number_of_perforations_;
// record the index of the first perforation
// TODO: it might not be needed if we refactor WellState to be a vector
// of states of individual well.
int first_perf_;
// well index for each perforation
std::vector<double> well_index_;
// depth for each perforation
std::vector<double> perf_depth_;
// reference depth for the BHP
double ref_depth_;
double well_efficiency_factor_;
// cell index for each well perforation
std::vector<int> well_cell_;
// saturation table nubmer for each well perforation
std::vector<int> saturation_table_number_;
const PhaseUsage* phase_usage_;
const std::vector<bool>* active_;
const VFPProperties* vfp_properties_;
double gravity_;
};
}
#include "WellInterface_impl.hpp"
#endif // OPM_WELLINTERFACE_HEADER_INCLUDED