Replaced TwoPhaseState -> SimulatorState

This commit is contained in:
Joakim Hove 2016-02-22 18:44:02 +01:00
parent 6ca10ea57d
commit df6e6a2347
2 changed files with 18 additions and 16 deletions

View File

@ -28,7 +28,7 @@
#include <opm/core/pressure/flow_bc.h>
#include <opm/core/linalg/LinearSolverInterface.hpp>
#include <opm/core/linalg/sparse_sys.h>
#include <opm/core/simulator/TwophaseState.hpp>
#include <opm/core/simulator/SimulatorState.hpp>
#include <opm/core/simulator/WellState.hpp>
#include <opm/common/ErrorMacros.hpp>
#include <opm/core/utility/miscUtilities.hpp>
@ -155,7 +155,7 @@ namespace Opm
/// May throw an exception if the number of iterations
/// exceed maxiter (set in constructor).
void IncompTpfa::solve(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state)
{
if (rock_comp_props_ != 0 && rock_comp_props_->isActive()) {
@ -169,7 +169,7 @@ namespace Opm
// Solve with no rock compressibility (linear eqn).
void IncompTpfa::solveIncomp(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state)
{
// Set up properties.
@ -207,7 +207,7 @@ namespace Opm
// Solve with rock compressibility (nonlinear eqn).
void IncompTpfa::solveRockComp(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state)
{
// This function is identical to CompressibleTpfa::solve().
@ -321,7 +321,7 @@ namespace Opm
/// Compute per-solve dynamic properties.
void IncompTpfa::computePerSolveDynamicData(const double /*dt*/,
const TwophaseState& state,
const SimulatorState& state,
const WellState& /*well_state*/)
{
// Computed here:
@ -369,7 +369,7 @@ namespace Opm
/// Compute per-iteration dynamic properties.
void IncompTpfa::computePerIterationDynamicData(const double /*dt*/,
const TwophaseState& state,
const SimulatorState& state,
const WellState& well_state)
{
// These are the variables that get computed by this function:
@ -396,7 +396,7 @@ namespace Opm
/// Compute the residual in h_->b and Jacobian in h_->A.
void IncompTpfa::assemble(const double dt,
const TwophaseState& state,
const SimulatorState& state,
const WellState& /*well_state*/)
{
const double* pressures = wells_ ? &pressures_[0] : &state.pressure()[0];
@ -462,7 +462,7 @@ namespace Opm
/// Compute the output.
void IncompTpfa::computeResults(TwophaseState& state,
void IncompTpfa::computeResults(SimulatorState& state,
WellState& well_state) const
{
// Make sure h_ contains the direct-solution matrix

View File

@ -20,6 +20,7 @@
#ifndef OPM_INCOMPTPFA_HEADER_INCLUDED
#define OPM_INCOMPTPFA_HEADER_INCLUDED
#include <opm/core/simulator/SimulatorState.hpp>
#include <opm/core/pressure/tpfa/ifs_tpfa.h>
#include <vector>
@ -34,8 +35,9 @@ namespace Opm
class IncompPropertiesInterface;
class RockCompressibility;
class LinearSolverInterface;
class TwophaseState;
class WellState;
class SimulatoreState;
/// Encapsulating a tpfa pressure solver for the incompressible-fluid case.
/// Supports gravity, wells controlled by bhp or reservoir rates,
@ -112,7 +114,7 @@ namespace Opm
/// May throw an exception if the number of iterations
/// exceed maxiter (set in constructor).
void solve(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state);
@ -122,28 +124,28 @@ namespace Opm
protected:
// Solve with no rock compressibility (linear eqn).
void solveIncomp(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state);
// Solve with rock compressibility (nonlinear eqn).
void solveRockComp(const double dt,
TwophaseState& state,
SimulatorState& state,
WellState& well_state);
private:
// Helper functions.
void computeStaticData();
virtual void computePerSolveDynamicData(const double dt,
const TwophaseState& state,
const SimulatorState& state,
const WellState& well_state);
void computePerIterationDynamicData(const double dt,
const TwophaseState& state,
const SimulatorState& state,
const WellState& well_state);
void assemble(const double dt,
const TwophaseState& state,
const SimulatorState& state,
const WellState& well_state);
void solveIncrement();
double residualNorm() const;
double incrementNorm() const;
void computeResults(TwophaseState& state,
void computeResults(SimulatorState& state,
WellState& well_state) const;
protected: