mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
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.
This commit is contained in:
@@ -321,7 +321,11 @@ namespace Opm {
|
|||||||
// Apply the update, with considering model-dependent limitations and
|
// Apply the update, with considering model-dependent limitations and
|
||||||
// chopping of the update.
|
// chopping of the update.
|
||||||
updateState(x,iteration);
|
updateState(x,iteration);
|
||||||
wellModel().updateWellState(xw, well_state);
|
|
||||||
|
if( nw > 0 )
|
||||||
|
{
|
||||||
|
wellModel().applySolutionWellState(x, well_state);
|
||||||
|
}
|
||||||
report.update_time += perfTimer.stop();
|
report.update_time += perfTimer.stop();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|||||||
@@ -141,6 +141,11 @@ namespace Opm
|
|||||||
// r = r - C D^-1 Rw
|
// r = r - C D^-1 Rw
|
||||||
virtual void apply(BVector& r) const;
|
virtual void apply(BVector& r) const;
|
||||||
|
|
||||||
|
// 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;
|
||||||
|
|
||||||
using WellInterface<TypeTag>::phaseUsage;
|
using WellInterface<TypeTag>::phaseUsage;
|
||||||
using WellInterface<TypeTag>::active;
|
using WellInterface<TypeTag>::active;
|
||||||
using WellInterface<TypeTag>::numberOfPerforations;
|
using WellInterface<TypeTag>::numberOfPerforations;
|
||||||
@@ -167,6 +172,9 @@ namespace Opm
|
|||||||
// TODO: maybe this function can go to some helper file.
|
// TODO: maybe this function can go to some helper file.
|
||||||
void localInvert(Mat& istlA) const;
|
void localInvert(Mat& istlA) const;
|
||||||
|
|
||||||
|
// xw = inv(D)*(rw - C*x)
|
||||||
|
void recoverSolutionWell(const BVector& x, BVector& xw) const;
|
||||||
|
|
||||||
// TODO: decide wether to use member function to refer to private member later
|
// TODO: decide wether to use member function to refer to private member later
|
||||||
using WellInterface<TypeTag>::vfp_properties_;
|
using WellInterface<TypeTag>::vfp_properties_;
|
||||||
using WellInterface<TypeTag>::gravity_;
|
using WellInterface<TypeTag>::gravity_;
|
||||||
|
|||||||
@@ -1904,4 +1904,36 @@ namespace Opm
|
|||||||
duneC_.mmtv(invDrw_, r);
|
duneC_.mmtv(invDrw_, r);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
template<typename TypeTag>
|
||||||
|
void
|
||||||
|
StandardWell<TypeTag>::
|
||||||
|
recoverSolutionWell(const BVector& x, BVector& xw) const
|
||||||
|
{
|
||||||
|
BVector resWell = resWell_;
|
||||||
|
// resWell = resWell - B * x
|
||||||
|
duneB_.mmv(x, resWell);
|
||||||
|
// xw = D^-1 * resWell
|
||||||
|
invDuneD_.mv(resWell, xw);
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
template<typename TypeTag>
|
||||||
|
void
|
||||||
|
StandardWell<TypeTag>::
|
||||||
|
applySolutionWellState(const BVector& x,
|
||||||
|
const ModelParameters& param,
|
||||||
|
WellState& well_state) const
|
||||||
|
{
|
||||||
|
BVector xw(1);
|
||||||
|
recoverSolutionWell(x, xw);
|
||||||
|
updateWellState(xw, param, well_state);
|
||||||
|
}
|
||||||
|
|
||||||
}
|
}
|
||||||
|
|||||||
@@ -98,6 +98,8 @@ enum WellVariablePositions {
|
|||||||
static const int solventSaturationIdx = BlackoilIndices::solventSaturationIdx;
|
static const int solventSaturationIdx = BlackoilIndices::solventSaturationIdx;
|
||||||
static const int polymerConcentrationIdx = BlackoilIndices::polymerConcentrationIdx;
|
static const int polymerConcentrationIdx = BlackoilIndices::polymerConcentrationIdx;
|
||||||
|
|
||||||
|
// TODO: where we should put these types, WellInterface or Well Model?
|
||||||
|
// or there is some other strategy, like TypeTag
|
||||||
typedef Dune::FieldVector<Scalar, numEq > VectorBlockType;
|
typedef Dune::FieldVector<Scalar, numEq > VectorBlockType;
|
||||||
typedef Dune::FieldMatrix<Scalar, numEq, numEq > MatrixBlockType;
|
typedef Dune::FieldMatrix<Scalar, numEq, numEq > MatrixBlockType;
|
||||||
typedef Dune::BCRSMatrix <MatrixBlockType> Mat;
|
typedef Dune::BCRSMatrix <MatrixBlockType> Mat;
|
||||||
@@ -175,8 +177,9 @@ enum WellVariablePositions {
|
|||||||
// apply well model with scaling of alpha
|
// apply well model with scaling of alpha
|
||||||
void applyScaleAdd(const Scalar alpha, const BVector& x, BVector& Ax) const;
|
void applyScaleAdd(const Scalar alpha, const BVector& x, BVector& Ax) const;
|
||||||
|
|
||||||
// xw = inv(D)*(rw - C*x)
|
// using the solution x to recover the solution xw for wells and applying
|
||||||
void recoverVariable(const BVector& x, BVector& xw) const;
|
// xw to update Well State
|
||||||
|
void applySolutionWellState(const BVector& x, WellState& well_state) const;
|
||||||
|
|
||||||
int flowPhaseToEbosCompIdx( const int phaseIdx ) const;
|
int flowPhaseToEbosCompIdx( const int phaseIdx ) const;
|
||||||
|
|
||||||
|
|||||||
@@ -431,6 +431,7 @@ namespace Opm {
|
|||||||
StandardWellsDense<TypeTag>::
|
StandardWellsDense<TypeTag>::
|
||||||
apply(const BVector& x, BVector& Ax) const
|
apply(const BVector& x, BVector& Ax) const
|
||||||
{
|
{
|
||||||
|
// TODO: do we still need localWellsActive()?
|
||||||
if ( ! localWellsActive() ) {
|
if ( ! localWellsActive() ) {
|
||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
@@ -486,28 +487,20 @@ namespace Opm {
|
|||||||
|
|
||||||
|
|
||||||
|
|
||||||
// xw = D^-1(resWell - B * x)
|
|
||||||
// TODO: this function should be moved to StandardWell
|
|
||||||
// xw should not appear in StandardWellsDense to avoid
|
|
||||||
// the data type to hold different types of xw
|
|
||||||
template<typename TypeTag>
|
template<typename TypeTag>
|
||||||
void
|
void
|
||||||
StandardWellsDense<TypeTag>::
|
StandardWellsDense<TypeTag>::
|
||||||
recoverVariable(const BVector& x, BVector& xw) const
|
applySolutionWellState(const BVector& x, WellState& well_state) const
|
||||||
{
|
{
|
||||||
if ( ! localWellsActive() ) {
|
for (auto& well : well_container_) {
|
||||||
return;
|
well->applySolutionWellState(x, param_, well_state);
|
||||||
}
|
}
|
||||||
BVector resWell = resWell_;
|
|
||||||
// resWell = resWell - B * x
|
|
||||||
duneB_.mmv(x, resWell);
|
|
||||||
// xw = D^-1 * resWell
|
|
||||||
invDuneD_.mv(resWell, xw);
|
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
template<typename TypeTag>
|
template<typename TypeTag>
|
||||||
int
|
int
|
||||||
StandardWellsDense<TypeTag>::
|
StandardWellsDense<TypeTag>::
|
||||||
|
|||||||
@@ -192,6 +192,11 @@ namespace Opm
|
|||||||
// r = r - C D^-1 Rw
|
// r = r - C D^-1 Rw
|
||||||
virtual void apply(BVector& r) const = 0;
|
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:
|
protected:
|
||||||
// TODO: some variables shared by all the wells should be made static
|
// TODO: some variables shared by all the wells should be made static
|
||||||
// well name
|
// well name
|
||||||
|
|||||||
Reference in New Issue
Block a user