2013-12-27 01:48:04 -06:00
|
|
|
/*
|
2014-02-24 19:52:10 -06:00
|
|
|
Copyright 2014 SINTEF ICT, Applied Mathematics.
|
2014-03-16 21:54:29 -05:00
|
|
|
Copyright 2014 STATOIL ASA.
|
2013-12-27 01:48:04 -06: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/>.
|
|
|
|
*/
|
|
|
|
|
2014-11-17 19:29:27 -06:00
|
|
|
#ifndef OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
|
|
|
|
#define OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
|
2013-12-27 01:48:04 -06:00
|
|
|
|
2014-03-14 01:12:26 -05:00
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
|
|
|
#include <opm/autodiff/AutoDiffHelpers.hpp>
|
2016-12-13 08:31:19 -06:00
|
|
|
#include <opm/autodiff/BlackoilModelEnums.hpp>
|
|
|
|
#include <opm/autodiff/BlackoilPropsAdFromDeck.hpp>
|
2014-09-26 01:54:57 -05:00
|
|
|
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
|
|
|
|
#include <opm/autodiff/LinearisedBlackoilResidual.hpp>
|
2013-12-27 01:48:04 -06:00
|
|
|
#include <opm/polymer/PolymerProperties.hpp>
|
2015-05-28 05:20:59 -05:00
|
|
|
#include <opm/polymer/fullyimplicit/WellStateFullyImplicitBlackoilPolymer.hpp>
|
2013-12-27 01:48:04 -06:00
|
|
|
#include <opm/polymer/fullyimplicit/PolymerPropsAd.hpp>
|
2015-06-22 02:43:01 -05:00
|
|
|
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
2017-02-10 09:07:25 -06:00
|
|
|
#include <opm/simulators/timestepping/SimulatorTimerInterface.hpp>
|
2017-03-16 04:22:27 -05:00
|
|
|
#include <opm/common/data/SimulationDataContainer.hpp>
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
struct UnstructuredGrid;
|
|
|
|
struct Wells;
|
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
class DerivedGeology;
|
|
|
|
class RockCompressibility;
|
2014-09-26 01:54:57 -05:00
|
|
|
class NewtonIterationBlackoilInterface;
|
2013-12-27 01:48:04 -06:00
|
|
|
class PolymerBlackoilState;
|
2014-10-28 21:25:06 -05:00
|
|
|
class WellStateFullyImplicitBlackoil;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
2014-02-24 19:52:10 -06:00
|
|
|
/// A fully implicit solver for the oil-water with polymer problem.
|
2013-12-27 01:48:04 -06:00
|
|
|
///
|
2014-02-24 19:52:10 -06:00
|
|
|
/// The simulator is capable of handling oil-water-polymer problems
|
|
|
|
/// It uses an industry-standard TPFA discretization with per-phase
|
2013-12-27 01:48:04 -06:00
|
|
|
/// upwind weighting of mobilities.
|
|
|
|
///
|
|
|
|
/// It uses automatic differentiation via the class AutoDiffBlock
|
|
|
|
/// to simplify assembly of the jacobian matrix.
|
|
|
|
class FullyImplicitCompressiblePolymerSolver
|
|
|
|
{
|
|
|
|
public:
|
2016-07-17 19:50:50 -05:00
|
|
|
typedef AutoDiffBlock<double> ADB;
|
|
|
|
typedef ADB::V V;
|
|
|
|
typedef ADB::M M;
|
|
|
|
typedef Eigen::Array<double,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::Dynamic,
|
|
|
|
Eigen::RowMajor> DataBlock;
|
2016-09-01 04:38:04 -05:00
|
|
|
|
2016-09-02 02:38:25 -05:00
|
|
|
struct ReservoirResidualQuant {
|
|
|
|
ReservoirResidualQuant();
|
|
|
|
std::vector<ADB> accum; // Accumulations
|
|
|
|
ADB mflux; // Mass flux (surface conditions)
|
|
|
|
ADB b; // Reciprocal FVF
|
|
|
|
ADB mu; // Viscosities
|
|
|
|
ADB rho; // Densities
|
|
|
|
ADB kr; // Permeabilities
|
|
|
|
ADB head; // Pressure drop across int. interfaces
|
|
|
|
ADB mob; // Phase mobility (per cell)
|
|
|
|
std::vector<ADB> ads; // Adsorption term.
|
|
|
|
};
|
|
|
|
|
2017-02-09 06:54:42 -06:00
|
|
|
struct SimulatorData : public Opm::FIPDataEnums {
|
2016-09-02 07:15:10 -05:00
|
|
|
SimulatorData(int num_phases)
|
|
|
|
: rq(num_phases)
|
2016-09-14 08:41:47 -05:00
|
|
|
, rsSat(ADB::null())
|
|
|
|
, rvSat(ADB::null())
|
2017-03-15 07:51:51 -05:00
|
|
|
, soMax()
|
2016-09-26 07:17:52 -05:00
|
|
|
, fip()
|
2016-09-02 07:15:10 -05:00
|
|
|
{
|
|
|
|
}
|
2016-09-26 07:17:52 -05:00
|
|
|
|
2017-02-09 06:54:42 -06:00
|
|
|
using Opm::FIPDataEnums::FipId;
|
|
|
|
using Opm::FIPDataEnums::fipValues;
|
|
|
|
|
2016-09-02 07:15:10 -05:00
|
|
|
std::vector<ReservoirResidualQuant> rq;
|
2016-09-14 08:41:47 -05:00
|
|
|
ADB rsSat;
|
|
|
|
ADB rvSat;
|
2017-03-15 07:51:51 -05:00
|
|
|
std::vector<double> soMax;
|
2017-02-09 06:54:42 -06:00
|
|
|
std::array<V, fipValues> fip;
|
2016-09-02 07:15:10 -05:00
|
|
|
};
|
2016-09-02 02:38:25 -05:00
|
|
|
|
2017-02-09 06:54:42 -06:00
|
|
|
typedef Opm::FIPData FIPDataType;
|
|
|
|
|
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
/// Construct a solver. It will retain references to the
|
|
|
|
/// arguments of this functions, and they are expected to
|
|
|
|
/// remain in scope for the lifetime of the solver.
|
|
|
|
/// \param[in] grid grid data structure
|
|
|
|
/// \param[in] fluid fluid properties
|
|
|
|
/// \param[in] geo rock properties
|
|
|
|
/// \param[in] rock_comp_props if non-null, rock compressibility properties
|
2014-02-24 19:52:10 -06:00
|
|
|
/// \param[in] polymer_props_ad polymer properties
|
2013-12-27 01:48:04 -06:00
|
|
|
/// \param[in] wells well structure
|
|
|
|
/// \param[in] linsolver linear solver
|
|
|
|
FullyImplicitCompressiblePolymerSolver(const UnstructuredGrid& grid ,
|
2016-12-29 09:35:24 -06:00
|
|
|
const BlackoilPropsAdFromDeck& fluid,
|
2016-06-28 01:38:20 -05:00
|
|
|
const DerivedGeology& geo ,
|
|
|
|
const RockCompressibility* rock_comp_props,
|
|
|
|
const PolymerPropsAd& polymer_props_ad,
|
|
|
|
const Wells& wells,
|
|
|
|
const NewtonIterationBlackoilInterface& linsolver);
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
/// Take a single forward step, modifiying
|
|
|
|
/// state.pressure()
|
|
|
|
/// state.faceflux()
|
|
|
|
/// state.saturation()
|
2014-02-24 19:52:10 -06:00
|
|
|
/// state.concentration()
|
2013-12-27 01:48:04 -06:00
|
|
|
/// wstate.bhp()
|
|
|
|
/// \param[in] dt time step size
|
|
|
|
/// \param[in] state reservoir state
|
|
|
|
/// \param[in] wstate well state
|
2014-02-24 19:52:10 -06:00
|
|
|
/// \param[in] polymer_inflow polymer influx
|
2015-05-28 05:20:59 -05:00
|
|
|
int
|
2016-07-05 05:20:19 -05:00
|
|
|
step(const SimulatorTimerInterface& timer,
|
2014-02-24 19:52:10 -06:00
|
|
|
PolymerBlackoilState& state ,
|
2015-05-28 05:20:59 -05:00
|
|
|
WellStateFullyImplicitBlackoilPolymer& wstate);
|
|
|
|
|
2016-07-13 05:06:54 -05:00
|
|
|
int linearizations() const;
|
2015-10-02 07:00:54 -05:00
|
|
|
int nonlinearIterations() const;
|
2015-05-28 06:15:42 -05:00
|
|
|
int linearIterations() const;
|
2016-06-20 19:36:15 -05:00
|
|
|
int wellIterations() const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
2015-06-10 06:22:43 -05:00
|
|
|
/// Not used by this class except to satisfy interface requirements.
|
|
|
|
typedef parameter::ParameterGroup SolverParameters;
|
|
|
|
|
2015-11-13 06:58:41 -06:00
|
|
|
/// There is no separate model class for this solver, return itself.
|
|
|
|
const FullyImplicitCompressiblePolymerSolver& model() const;
|
|
|
|
|
|
|
|
/// Evaluate the relative changes in the physical variables.
|
|
|
|
double relativeChange(const PolymerBlackoilState& previous,
|
|
|
|
const PolymerBlackoilState& current ) const;
|
|
|
|
|
2016-09-02 07:15:10 -05:00
|
|
|
/// Return reservoir simulation data (for output functionality)
|
2017-03-16 04:22:27 -05:00
|
|
|
const SimulatorData& getSimulatorData(const SimulationDataContainer&) const {
|
2016-09-02 07:15:10 -05:00
|
|
|
return sd_;
|
2016-09-02 02:38:25 -05:00
|
|
|
}
|
|
|
|
|
2017-02-09 06:54:42 -06:00
|
|
|
/// Return reservoir simulation data (for output functionality)
|
|
|
|
FIPDataType getFIPData() const {
|
|
|
|
return FIPDataType( sd_.fip );
|
|
|
|
}
|
|
|
|
|
2016-07-17 19:50:50 -05:00
|
|
|
/// Compute fluid in place.
|
2016-08-09 02:30:25 -05:00
|
|
|
/// \param[in] ReservoirState
|
|
|
|
/// \param[in] WellState
|
|
|
|
/// \param[in] FIPNUM for active cells not global cells.
|
2017-02-09 06:54:42 -06:00
|
|
|
/// \return fluid in place, number of fip regions, each region contains 5 values which are liquid, vapour, water, free gas and dissolved gas.
|
|
|
|
std::vector<std::vector<double> >
|
2016-07-18 20:33:03 -05:00
|
|
|
computeFluidInPlace(const PolymerBlackoilState& x,
|
|
|
|
const std::vector<int>& fipnum);
|
2016-07-17 19:50:50 -05:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
private:
|
2016-07-17 19:50:50 -05:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
struct SolutionState {
|
|
|
|
SolutionState(const int np);
|
|
|
|
ADB pressure;
|
2014-11-20 11:57:29 -06:00
|
|
|
ADB temperature;
|
2013-12-27 01:48:04 -06:00
|
|
|
std::vector<ADB> saturation;
|
|
|
|
ADB concentration;
|
|
|
|
ADB qs;
|
|
|
|
ADB bhp;
|
|
|
|
};
|
|
|
|
|
|
|
|
struct WellOps {
|
|
|
|
WellOps(const Wells& wells);
|
2015-09-08 02:57:17 -05:00
|
|
|
Eigen::SparseMatrix<double> w2p; // well -> perf (scatter)
|
|
|
|
Eigen::SparseMatrix<double> p2w; // perf -> well (gather)
|
2013-12-27 01:48:04 -06:00
|
|
|
};
|
|
|
|
|
2016-12-13 08:31:19 -06:00
|
|
|
enum { Water = Opm::Water,
|
|
|
|
Oil = Opm::Oil };
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
// Member data
|
|
|
|
const UnstructuredGrid& grid_;
|
2016-12-29 09:35:24 -06:00
|
|
|
const BlackoilPropsAdFromDeck& fluid_;
|
2013-12-27 01:48:04 -06:00
|
|
|
const DerivedGeology& geo_;
|
|
|
|
const RockCompressibility* rock_comp_props_;
|
|
|
|
const PolymerPropsAd& polymer_props_ad_;
|
|
|
|
const Wells& wells_;
|
2014-09-26 01:54:57 -05:00
|
|
|
const NewtonIterationBlackoilInterface& linsolver_;
|
2013-12-27 01:48:04 -06:00
|
|
|
const std::vector<int> cells_; // All grid cells
|
|
|
|
HelperOps ops_;
|
|
|
|
const WellOps wops_;
|
|
|
|
const M grav_;
|
2014-01-20 02:25:30 -06:00
|
|
|
V cmax_;
|
2014-09-24 22:16:56 -05:00
|
|
|
std::vector<PhasePresence> phaseCondition_;
|
2016-09-02 07:15:10 -05:00
|
|
|
SimulatorData sd_;
|
2013-12-27 01:48:04 -06:00
|
|
|
// The mass_balance vector has one element for each active phase,
|
|
|
|
// each of which has size equal to the number of cells.
|
|
|
|
// The well_eq has size equal to the number of wells.
|
2014-09-26 01:54:57 -05:00
|
|
|
LinearisedBlackoilResidual residual_;
|
2015-05-28 05:20:59 -05:00
|
|
|
|
2016-07-13 05:06:54 -05:00
|
|
|
unsigned int linearizations_;
|
2015-05-28 05:20:59 -05:00
|
|
|
unsigned int newtonIterations_;
|
|
|
|
unsigned int linearIterations_;
|
2016-06-20 19:36:15 -05:00
|
|
|
unsigned int wellIterations_;
|
2015-05-28 05:20:59 -05:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
// Private methods.
|
|
|
|
SolutionState
|
|
|
|
constantState(const PolymerBlackoilState& x,
|
2014-10-28 21:25:06 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& xw);
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
SolutionState
|
|
|
|
variableState(const PolymerBlackoilState& x,
|
2014-10-28 21:25:06 -05:00
|
|
|
const WellStateFullyImplicitBlackoil& xw);
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
void
|
|
|
|
computeAccum(const SolutionState& state,
|
|
|
|
const int aix );
|
|
|
|
|
|
|
|
void
|
2014-02-24 19:52:10 -06:00
|
|
|
assemble(const double dt,
|
2013-12-27 01:48:04 -06:00
|
|
|
const PolymerBlackoilState& x,
|
2017-02-09 06:54:42 -06:00
|
|
|
const WellStateFullyImplicitBlackoil& xw,
|
2014-10-29 00:38:15 -05:00
|
|
|
const std::vector<double>& polymer_inflow);
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
V solveJacobianSystem() const;
|
|
|
|
|
|
|
|
void updateState(const V& dx,
|
|
|
|
PolymerBlackoilState& state,
|
2014-10-28 21:25:06 -05:00
|
|
|
WellStateFullyImplicitBlackoil& well_state) const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
std::vector<ADB>
|
|
|
|
computeRelPerm(const SolutionState& state) const;
|
|
|
|
|
|
|
|
std::vector<ADB>
|
2015-05-29 09:31:32 -05:00
|
|
|
computePressures(const SolutionState& state) const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
void
|
|
|
|
computeMassFlux(const int actph ,
|
|
|
|
const V& transi,
|
|
|
|
const std::vector<ADB>& kr ,
|
|
|
|
const SolutionState& state );
|
2014-02-24 19:52:10 -06:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
void
|
|
|
|
computeMassFlux(const V& trans,
|
|
|
|
const ADB& mc,
|
|
|
|
const ADB& kro,
|
|
|
|
const ADB& krw_eff,
|
|
|
|
const SolutionState& state);
|
|
|
|
|
|
|
|
std::vector<ADB>
|
|
|
|
computeFracFlow(const ADB& kro,
|
|
|
|
const ADB& krw_eff,
|
|
|
|
const ADB& c) const;
|
2014-02-24 19:52:10 -06:00
|
|
|
|
2014-01-20 02:25:30 -06:00
|
|
|
void
|
2014-12-07 23:45:30 -06:00
|
|
|
computeCmax(PolymerBlackoilState& state);
|
2014-02-24 19:52:10 -06:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
ADB
|
|
|
|
computeMc(const SolutionState& state) const;
|
2014-02-24 19:52:10 -06:00
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
ADB
|
2014-02-24 19:52:10 -06:00
|
|
|
rockPorosity(const ADB& p) const;
|
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
ADB
|
2014-02-24 19:52:10 -06:00
|
|
|
rockPermeability(const ADB& p) const;
|
|
|
|
|
2013-12-27 01:48:04 -06:00
|
|
|
double
|
|
|
|
residualNorm() const;
|
|
|
|
|
|
|
|
ADB
|
2014-09-24 22:16:56 -05:00
|
|
|
fluidViscosity(const int phase,
|
|
|
|
const ADB& p ,
|
2014-11-20 11:57:29 -06:00
|
|
|
const ADB& T ,
|
2014-09-24 22:16:56 -05:00
|
|
|
const std::vector<PhasePresence>& cond,
|
|
|
|
const std::vector<int>& cells) const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
ADB
|
2014-09-24 22:16:56 -05:00
|
|
|
fluidReciprocFVF(const int phase,
|
|
|
|
const ADB& p ,
|
2014-11-20 11:57:29 -06:00
|
|
|
const ADB& T ,
|
2014-09-26 01:06:01 -05:00
|
|
|
const std::vector<PhasePresence>& cond,
|
2014-09-24 22:16:56 -05:00
|
|
|
const std::vector<int>& cells) const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
ADB
|
2014-09-24 22:16:56 -05:00
|
|
|
fluidDensity(const int phase,
|
|
|
|
const ADB& p ,
|
2014-11-20 11:57:29 -06:00
|
|
|
const ADB& T ,
|
2014-09-24 22:16:56 -05:00
|
|
|
const std::vector<PhasePresence>& cond,
|
|
|
|
const std::vector<int>& cells) const;
|
2013-12-27 01:48:04 -06:00
|
|
|
|
|
|
|
ADB
|
|
|
|
poroMult(const ADB& p) const;
|
|
|
|
|
|
|
|
ADB
|
|
|
|
transMult(const ADB& p) const;
|
2017-02-09 06:54:42 -06:00
|
|
|
|
2014-09-24 22:16:56 -05:00
|
|
|
const std::vector<PhasePresence>
|
2014-09-26 01:06:01 -05:00
|
|
|
phaseCondition() const { return phaseCondition_; }
|
2017-02-09 06:54:42 -06:00
|
|
|
|
2014-09-24 22:16:56 -05:00
|
|
|
void
|
|
|
|
classifyCondition(const PolymerBlackoilState& state);
|
2013-12-27 01:48:04 -06:00
|
|
|
};
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
|
2014-11-17 19:29:27 -06:00
|
|
|
#endif // OPM_FULLYIMPLICITCOMPRESSIBLEPOLYMERSOLVER_HEADER_INCLUDED
|