ImplicitTransportDefs.hpp
Go to the documentation of this file.
1#ifndef OPENRS_IMPLICITTRANSPORTDEFS_HEADER
2#define OPENRS_IMPLICITTRANSPORTDEFS_HEADER
3
4#include <opm/common/utility/platform_dependent/disable_warnings.h>
5
6#include <dune/common/fvector.hh>
7#include <dune/common/fmatrix.hh>
8
9#include <dune/istl/operators.hh>
10#include <dune/istl/solvers.hh>
11#include <dune/istl/bvector.hh>
12#include <dune/istl/bcrsmatrix.hh>
13#include <dune/istl/preconditioners.hh>
14
15#include <opm/common/utility/platform_dependent/reenable_warnings.h>
16
17#include <opm/grid/common/GridAdapter.hpp>
18
19#include <opm/grid/UnstructuredGrid.h>
22#include <opm/common/utility/parameters/ParameterGroup.hpp>
23
25
26#include <algorithm>
27#include <cstddef>
28#include <memory>
29#include <string>
30#include <vector>
31
32namespace Opm {
33 template <class Vector>
34 class MaxNormStl {
35 public:
36 static double
37 norm(const Vector& v) {
39 }
40 };
41
42 template <class Vector>
44 public:
45 static double
46 norm(const Vector& v) {
47 return v.infinity_norm();
48 }
49 };
50
51 template <int np = 2>
53 public:
54 ReservoirState(const UnstructuredGrid* g)
55 : press_ (g->number_of_cells),
56 fpress_(g->number_of_faces),
57 flux_ (g->number_of_faces),
58 sat_ (np * g->number_of_cells)
59 {}
60
61 ::std::vector<double>& pressure () { return press_ ; }
62 ::std::vector<double>& facepressure() { return fpress_; }
63 ::std::vector<double>& faceflux () { return flux_ ; }
64 ::std::vector<double>& saturation () { return sat_ ; }
65
66 const ::std::vector<double>& faceflux () const { return flux_; }
67 const ::std::vector<double>& saturation () const { return sat_ ; }
68
69 private:
70 ::std::vector<double> press_ ;
71 ::std::vector<double> fpress_;
72 ::std::vector<double> flux_ ;
73 ::std::vector<double> sat_ ;
74 };
75
76 class Rock {
77 public:
78 Rock(::std::size_t nc, ::std::size_t dim)
79 : dim_ (dim ),
80 perm_(nc * dim * dim),
81 poro_(nc ) {}
82
83 const ::std::vector<double>& perm() const { return perm_; }
84 const ::std::vector<double>& poro() const { return poro_; }
85
86 void
87 perm_homogeneous(double k) {
88 setVector(0.0, perm_);
89
90 const ::std::size_t d2 = dim_ * dim_;
91
92 for (::std::size_t c = 0, nc = poro_.size(); c < nc; ++c) {
93 for (::std::size_t i = 0; i < dim_; ++i) {
94 perm_[c*d2 + i*(dim_ + 1)] = k;
95 }
96 }
97 }
98
99 void
100 poro_homogeneous(double phi) {
101 setVector(phi, poro_);
102 }
103
104 private:
105 void
106 setVector(double x, ::std::vector<double>& v) {
107 ::std::fill(v.begin(), v.end(), x);
108 }
109
110 ::std::size_t dim_ ;
111 ::std::vector<double> perm_;
112 ::std::vector<double> poro_;
113 };
114
116 public:
118 : r_(r)
119 {}
120
121 template <class Sat ,
122 class Mob ,
123 class DMob>
124 void
125 mobility(int c, const Sat& s, Mob& mob, DMob& dmob) const {
126 const double s1 = s[0];
127
128 r_.phaseMobilities (c, s1, mob );
129 r_.phaseMobilitiesDeriv(c, s1, dmob);
130 }
131
132 template <class Sat ,
133 class PC ,
134 class DPC>
135 void
136 pc(int c, const Sat& s, PC& pc_arg, DPC& dpc) const {
137 const double s1 = s[0];
138 pc_arg = r_.capillaryPressure(c, s1);
139 dpc = r_.capillaryPressureDeriv(c, s1);
140 }
141
142 double density(int p) const {
143 if (p == 0) {
144 return r_.densityFirstPhase();
145 } else {
146 return r_.densitySecondPhase();
147 }
148 }
149
150 double s_min(int c) const {
151 return r_.s_min(c);
152 }
153
154 double s_max(int c) const {
155 return r_.s_max(c);
156 }
157
158 private:
160 };
161
163 public:
164 typedef Dune::FieldVector<double, 1> ScalarVectorBlockType;
165 typedef Dune::FieldMatrix<double, 1, 1> ScalarMatrixBlockType;
166
167 typedef Dune::BlockVector<ScalarVectorBlockType> ScalarBlockVector;
168 typedef Dune::BCRSMatrix <ScalarMatrixBlockType> ScalarBCRSMatrix;
169
170 void
172 const ScalarBlockVector& b,
174 {
175 Dune::MatrixAdapter<ScalarBCRSMatrix,
177 ScalarBlockVector> opA(A);
178
179 Dune::SeqILU<ScalarBCRSMatrix,
181 ScalarBlockVector> precond(A, 1.0);
182
183 int maxit = A.N();
184 double tol = 5.0e-7;
185 int verb = 0;
186
187 Dune::BiCGSTABSolver<ScalarBlockVector>
188 solver(opA, precond, tol, maxit, verb);
189
190 ScalarBlockVector bcpy(b);
191 Dune::InverseOperatorResult res;
192 solver.apply(x, bcpy, res);
193 }
194 };
195
197 public:
198 TransportSource() : nsrc(0), pf(0) {}
199
200 int nsrc;
201 int pf;
202 ::std::vector< int > cell;
203 ::std::vector<double> pressure;
204 ::std::vector<double> flux;
205 ::std::vector<double> saturation;
206 ::std::vector<double> periodic_cells;
207 ::std::vector<double> periodic_faces;
208 };
209
210 template <class Arr>
211 void
212 append_transport_source(int c, double p, double v, const Arr& s,
213 TransportSource& src)
214 {
215 src.cell .push_back(c);
216 src.pressure .push_back(p);
217 src.flux .push_back(v);
218 src.saturation.insert(src.saturation.end(),
219 s.begin(), s.end());
220 ++src.nsrc;
221 }
222
223 void
224 compute_porevolume(const UnstructuredGrid* g,
225 const Rock& rock,
226 std::vector<double>& porevol);
227} // namespace Opm
228
229#endif // OPENRS_IMPLICITTRANSPORTDEFS_HEADER
Definition: ImplicitTransportDefs.hpp:162
Dune::BCRSMatrix< ScalarMatrixBlockType > ScalarBCRSMatrix
Definition: ImplicitTransportDefs.hpp:168
void solve(const ScalarBCRSMatrix &A, const ScalarBlockVector &b, ScalarBlockVector &x)
Definition: ImplicitTransportDefs.hpp:171
Dune::FieldMatrix< double, 1, 1 > ScalarMatrixBlockType
Definition: ImplicitTransportDefs.hpp:165
Dune::FieldVector< double, 1 > ScalarVectorBlockType
Definition: ImplicitTransportDefs.hpp:164
Dune::BlockVector< ScalarVectorBlockType > ScalarBlockVector
Definition: ImplicitTransportDefs.hpp:167
Definition: ImplicitTransportDefs.hpp:43
static double norm(const Vector &v)
Definition: ImplicitTransportDefs.hpp:46
Definition: ImplicitTransportDefs.hpp:34
static double norm(const Vector &v)
Definition: ImplicitTransportDefs.hpp:37
void phaseMobilities(int cell_index, double saturation, Vector &mobility) const
Mobilities for both phases.
Definition: ReservoirPropertyCapillary_impl.hpp:93
void phaseMobilitiesDeriv(int c, double s, Vector &dmob) const
Definition: ReservoirPropertyCapillary_impl.hpp:104
double densitySecondPhase() const
Density of second (oil) phase.
Definition: ReservoirPropertyCommon_impl.hpp:373
double densityFirstPhase() const
Density of first (water) phase.
Definition: ReservoirPropertyCommon_impl.hpp:366
double capillaryPressure(int cell_index, double saturation) const
Capillary pressure.
Definition: ReservoirPropertyCommon_impl.hpp:467
double capillaryPressureDeriv(int c, double s) const
Derivative of Capillary pressure.
Definition: ReservoirPropertyCommon_impl.hpp:480
double s_min(int c) const
Definition: ReservoirPropertyCommon_impl.hpp:493
double s_max(int c) const
Definition: ReservoirPropertyCommon_impl.hpp:506
Definition: ImplicitTransportDefs.hpp:52
ReservoirState(const UnstructuredGrid *g)
Definition: ImplicitTransportDefs.hpp:54
::std::vector< double > & pressure()
Definition: ImplicitTransportDefs.hpp:61
::std::vector< double > & faceflux()
Definition: ImplicitTransportDefs.hpp:63
::std::vector< double > & saturation()
Definition: ImplicitTransportDefs.hpp:64
const ::std::vector< double > & saturation() const
Definition: ImplicitTransportDefs.hpp:67
const ::std::vector< double > & faceflux() const
Definition: ImplicitTransportDefs.hpp:66
::std::vector< double > & facepressure()
Definition: ImplicitTransportDefs.hpp:62
A property class for porous media rock.
Definition: ImplicitTransportDefs.hpp:76
void poro_homogeneous(double phi)
Definition: ImplicitTransportDefs.hpp:100
void perm_homogeneous(double k)
Definition: ImplicitTransportDefs.hpp:87
Rock(::std::size_t nc, ::std::size_t dim)
Definition: ImplicitTransportDefs.hpp:78
const ::std::vector< double > & poro() const
Definition: ImplicitTransportDefs.hpp:84
const ::std::vector< double > & perm() const
Definition: ImplicitTransportDefs.hpp:83
Definition: ImplicitTransportDefs.hpp:196
::std::vector< double > periodic_faces
Definition: ImplicitTransportDefs.hpp:207
int pf
Definition: ImplicitTransportDefs.hpp:201
::std::vector< double > saturation
Definition: ImplicitTransportDefs.hpp:205
::std::vector< int > cell
Definition: ImplicitTransportDefs.hpp:202
TransportSource()
Definition: ImplicitTransportDefs.hpp:198
::std::vector< double > periodic_cells
Definition: ImplicitTransportDefs.hpp:206
int nsrc
Definition: ImplicitTransportDefs.hpp:200
::std::vector< double > flux
Definition: ImplicitTransportDefs.hpp:204
::std::vector< double > pressure
Definition: ImplicitTransportDefs.hpp:203
Definition: ImplicitTransportDefs.hpp:115
double s_max(int c) const
Definition: ImplicitTransportDefs.hpp:154
TwophaseFluidWrapper(const Opm::ReservoirPropertyCapillary< 3 > &r)
Definition: ImplicitTransportDefs.hpp:117
double s_min(int c) const
Definition: ImplicitTransportDefs.hpp:150
double density(int p) const
Definition: ImplicitTransportDefs.hpp:142
void pc(int c, const Sat &s, PC &pc_arg, DPC &dpc) const
Definition: ImplicitTransportDefs.hpp:136
void mobility(int c, const Sat &s, Mob &mob, DMob &dmob) const
Definition: ImplicitTransportDefs.hpp:125
Dune::BlockVector< Dune::FieldVector< double, 1 > > Vector
A vector holding our RHS.
Definition: matrixops.hpp:33
Definition: ImplicitAssembly.hpp:43
void compute_porevolume(const UnstructuredGrid *g, const Rock &rock, std::vector< double > &porevol)
void append_transport_source(int c, double p, double v, const Arr &s, TransportSource &src)
Definition: ImplicitTransportDefs.hpp:212