CompressibleTpfa.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2012 SINTEF ICT, Applied Mathematics.
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef OPM_COMPRESSIBLETPFA_HEADER_INCLUDED
21 #define OPM_COMPRESSIBLETPFA_HEADER_INCLUDED
22 
23 
24 #include <vector>
25 
26 struct UnstructuredGrid;
27 struct cfs_tpfa_res_data;
28 struct Wells;
30 
31 namespace Opm
32 {
33 
34  class BlackoilState;
35  class BlackoilPropertiesInterface;
36  class RockCompressibility;
37  class LinearSolverInterface;
38  class WellState;
39 
45  {
46  public:
65  const BlackoilPropertiesInterface& props,
66  const RockCompressibility* rock_comp_props,
67  const LinearSolverInterface& linsolver,
68  const double residual_tol,
69  const double change_tol,
70  const int maxiter,
71  const double* gravity,
72  const Wells* wells);
73 
75  virtual ~CompressibleTpfa();
76 
80  void solve(const double dt,
81  BlackoilState& state,
82  WellState& well_state);
83 
90  bool singularPressure() const;
91 
92  private:
93  virtual void computePerSolveDynamicData(const double dt,
94  const BlackoilState& state,
95  const WellState& well_state);
96  void computePerIterationDynamicData(const double dt,
97  const BlackoilState& state,
98  const WellState& well_state);
99  virtual void computeCellDynamicData(const double dt,
100  const BlackoilState& state,
101  const WellState& well_state);
102  void computeFaceDynamicData(const double dt,
103  const BlackoilState& state,
104  const WellState& well_state);
105  void computeWellDynamicData(const double dt,
106  const BlackoilState& state,
107  const WellState& well_state);
108  void assemble(const double dt,
109  const BlackoilState& state,
110  const WellState& well_state);
111  void solveIncrement();
112  double residualNorm() const;
113  double incrementNorm() const;
114  void computeResults(BlackoilState& state,
115  WellState& well_state) const;
116  protected:
117  void computeWellPotentials(const BlackoilState& state);
118 
119  // ------ Data that will remain unmodified after construction. ------
124  const double residual_tol_;
125  const double change_tol_;
126  const int maxiter_;
127  const double* gravity_; // May be NULL
128  const Wells* wells_; // May be NULL, outside may modify controls (only) between calls to solve().
129  std::vector<double> htrans_;
130  std::vector<double> trans_ ;
131  std::vector<int> allcells_;
132 
133  // ------ Internal data for the cfs_tpfa_res solver. ------
135 
136  // ------ Data that will be modified for every solve. ------
137  std::vector<double> wellperf_wdp_;
138  std::vector<double> initial_porevol_;
139 
140  // ------ Data that will be modified for every solver iteration. ------
141  std::vector<double> cell_A_;
142  std::vector<double> cell_dA_;
143  std::vector<double> cell_viscosity_;
144  std::vector<double> cell_phasemob_;
145  std::vector<double> cell_voldisc_;
146  std::vector<double> face_A_;
147  std::vector<double> face_phasemob_;
148  std::vector<double> face_gravcap_;
149  std::vector<double> wellperf_A_;
150  std::vector<double> wellperf_phasemob_;
151  std::vector<double> porevol_; // Only modified if rock_comp_props_ is non-null.
152  std::vector<double> rock_comp_; // Empty unless rock_comp_props_ is non-null.
153  // The update to be applied to the pressures (cell and bhp).
154  std::vector<double> pressure_increment_;
155  // True if the matrix assembled would be singular but for the
156  // adjustment made in the cfs_*_assemble() calls. This happens
157  // if everything is incompressible and there are no pressure
158  // conditions.
159  bool singular_;
160  };
161 
162 } // namespace Opm
163 
164 
165 #endif // OPM_COMPRESSIBLETPFA_HEADER_INCLUDED
Definition: wells.h:50
struct cfs_tpfa_res_data * h_
Definition: CompressibleTpfa.hpp:134
std::vector< double > wellperf_phasemob_
Definition: CompressibleTpfa.hpp:150
Definition: grid.h:98
std::vector< double > face_A_
Definition: CompressibleTpfa.hpp:146
bool singularPressure() const
After solve(), was the resulting pressure singular. Returns true if the pressure is singular in the f...
Definition: flow_bc.h:39
std::vector< double > face_gravcap_
Definition: CompressibleTpfa.hpp:148
Definition: AnisotropicEikonal.hpp:43
std::vector< double > cell_dA_
Definition: CompressibleTpfa.hpp:142
const double change_tol_
Definition: CompressibleTpfa.hpp:125
Definition: cfs_tpfa_residual.h:98
std::vector< double > cell_viscosity_
Definition: CompressibleTpfa.hpp:143
std::vector< double > trans_
Definition: CompressibleTpfa.hpp:130
const double * gravity_
Definition: CompressibleTpfa.hpp:127
std::vector< double > cell_phasemob_
Definition: CompressibleTpfa.hpp:144
std::vector< double > initial_porevol_
Definition: CompressibleTpfa.hpp:138
void solve(const double dt, BlackoilState &state, WellState &well_state)
const BlackoilPropertiesInterface & props_
Definition: CompressibleTpfa.hpp:121
std::vector< double > htrans_
Definition: CompressibleTpfa.hpp:129
Definition: BlackoilPropertiesInterface.hpp:37
bool singular_
Definition: CompressibleTpfa.hpp:159
std::vector< double > face_phasemob_
Definition: CompressibleTpfa.hpp:147
Definition: RockCompressibility.hpp:33
void computeWellPotentials(const BlackoilState &state)
const Wells * wells_
Definition: CompressibleTpfa.hpp:128
Simulator state for a blackoil simulator.
Definition: BlackoilState.hpp:33
const LinearSolverInterface & linsolver_
Definition: CompressibleTpfa.hpp:123
The state of a set of wells.
Definition: WellState.hpp:36
std::vector< double > porevol_
Definition: CompressibleTpfa.hpp:151
std::vector< double > wellperf_A_
Definition: CompressibleTpfa.hpp:149
std::vector< double > wellperf_wdp_
Definition: CompressibleTpfa.hpp:137
const double gravity
Definition: Units.hpp:120
const double residual_tol_
Definition: CompressibleTpfa.hpp:124
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
const UnstructuredGrid & grid_
Definition: CompressibleTpfa.hpp:120
virtual ~CompressibleTpfa()
Destructor.
std::vector< double > rock_comp_
Definition: CompressibleTpfa.hpp:152
CompressibleTpfa(const UnstructuredGrid &grid, const BlackoilPropertiesInterface &props, const RockCompressibility *rock_comp_props, const LinearSolverInterface &linsolver, const double residual_tol, const double change_tol, const int maxiter, const double *gravity, const Wells *wells)
Definition: CompressibleTpfa.hpp:44
std::vector< double > pressure_increment_
Definition: CompressibleTpfa.hpp:154
Abstract interface for linear solvers.
Definition: LinearSolverInterface.hpp:32
std::vector< int > allcells_
Definition: CompressibleTpfa.hpp:131
const int maxiter_
Definition: CompressibleTpfa.hpp:126
std::vector< double > cell_A_
Definition: CompressibleTpfa.hpp:141
std::vector< double > cell_voldisc_
Definition: CompressibleTpfa.hpp:145
const RockCompressibility * rock_comp_props_
Definition: CompressibleTpfa.hpp:122