BlackoilPolymerModel.hpp
Go to the documentation of this file.
1/*
2 Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
3 Copyright 2014 STATOIL ASA.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_BLACKOILPOLYMERMODEL_HEADER_INCLUDED
22#define OPM_BLACKOILPOLYMERMODEL_HEADER_INCLUDED
23
24#include <opm/autodiff/BlackoilModelBase.hpp>
25#include <opm/autodiff/BlackoilModelParameters.hpp>
30
31namespace Opm {
32
42 template<class Grid>
43 class BlackoilPolymerModel : public BlackoilModelBase<Grid, BlackoilPolymerModel<Grid> >
44 {
45 public:
46
47 // --------- Types and enums ---------
48
49 typedef BlackoilModelBase<Grid, BlackoilPolymerModel<Grid> > Base;
50 typedef typename Base::ReservoirState ReservoirState;
51 typedef typename Base::WellState WellState;
52 // The next line requires C++11 support available in g++ 4.7.
53 // friend Base;
54 friend class BlackoilModelBase<Grid, BlackoilPolymerModel<Grid> >;
55
75 BlackoilPolymerModel(const typename Base::ModelParameters& param,
76 const Grid& grid,
77 const BlackoilPropsAdInterface& fluid,
78 const DerivedGeology& geo,
79 const RockCompressibility* rock_comp_props,
80 const PolymerPropsAd& polymer_props_ad,
81 const Wells* wells,
82 const NewtonIterationBlackoilInterface& linsolver,
83 EclipseStateConstPtr eclipse_state,
84 const bool has_disgas,
85 const bool has_vapoil,
86 const bool has_polymer,
87 const bool has_plyshlog,
88 const bool has_shrate,
89 const std::vector<double>& wells_rep_radius,
90 const std::vector<double>& wells_perf_length,
91 const std::vector<double>& wells_bore_diameter,
92 const bool terminal_output);
93
98 void prepareStep(const double dt,
99 ReservoirState& reservoir_state,
100 WellState& well_state);
101
106 void afterStep(const double dt,
107 ReservoirState& reservoir_state,
108 WellState& well_state);
109
114 void updateState(const V& dx,
115 ReservoirState& reservoir_state,
116 WellState& well_state);
117
122 void assemble(const ReservoirState& reservoir_state,
123 WellState& well_state,
124 const bool initial_assembly);
125
126
127 protected:
128
129 // --------- Types and enums ---------
130
131 typedef typename Base::SolutionState SolutionState;
132 typedef typename Base::DataBlock DataBlock;
133 enum { Concentration = CanonicalVariablePositions::Next };
134
135 // --------- Data members ---------
136
138 const bool has_polymer_;
139 const bool has_plyshlog_;
140 const bool has_shrate_;
141 const int poly_pos_;
143
144 // representative radius and perforation length of well perforations
145 // to be used in shear-thinning computation.
146 std::vector<double> wells_rep_radius_;
147 std::vector<double> wells_perf_length_;
148 // wellbore diameters
149 std::vector<double> wells_bore_diameter_;
150
151 // shear-thinning factor for cell faces
152 std::vector<double> shear_mult_faces_;
153 // shear-thinning factor for well perforations
154 std::vector<double> shear_mult_wells_;
155
156 // Need to declare Base members we want to use here.
157 using Base::grid_;
158 using Base::fluid_;
159 using Base::geo_;
160 using Base::rock_comp_props_;
161 using Base::wells_;
162 using Base::linsolver_;
163 using Base::active_;
164 using Base::canph_;
165 using Base::cells_;
166 using Base::ops_;
167 using Base::wops_;
168 using Base::has_disgas_;
169 using Base::has_vapoil_;
170 using Base::param_;
171 using Base::use_threshold_pressure_;
172 using Base::threshold_pressures_by_interior_face_;
173 using Base::rq_;
174 using Base::phaseCondition_;
175 using Base::well_perforation_pressure_diffs_;
176 using Base::residual_;
177 using Base::terminal_output_;
178 using Base::primalVariable_;
179 using Base::pvdt_;
180
181 // --------- Protected methods ---------
182
183 // Need to declare Base members we want to use here.
184 using Base::wellsActive;
185 using Base::wells;
186 using Base::variableState;
187 using Base::computePressures;
188 using Base::computeGasPressure;
189 using Base::applyThresholdPressures;
190 using Base::fluidViscosity;
191 using Base::fluidReciprocFVF;
192 using Base::fluidDensity;
193 using Base::fluidRsSat;
194 using Base::fluidRvSat;
195 using Base::poroMult;
196 using Base::transMult;
197 using Base::updatePrimalVariableFromState;
198 using Base::updatePhaseCondFromPrimalVariable;
199 using Base::dpMaxRel;
200 using Base::dsMax;
201 using Base::drMaxRel;
202 using Base::maxResidualAllowed;
203
204 using Base::updateWellControls;
205 using Base::computeWellConnectionPressures;
206 using Base::addWellControlEq;
207 using Base::computeRelPerm;
208
209
210 void
211 makeConstantState(SolutionState& state) const;
212
213 std::vector<V>
215 const WellState& xw) const;
216
217 std::vector<int>
218 variableStateIndices() const;
219
222 const std::vector<int>& indices,
223 std::vector<ADB>& vars) const;
224
225 void
226 computeAccum(const SolutionState& state,
227 const int aix );
228
229 void
231
232 void
233 addWellContributionToMassBalanceEq(const std::vector<ADB>& cq_s,
234 const SolutionState& state,
235 WellState& xw);
236
237 void
238 computeMassFlux(const int actph ,
239 const V& transi,
240 const ADB& kr ,
241 const ADB& p ,
242 const SolutionState& state );
243
244 void
246
247 ADB
248 computeMc(const SolutionState& state) const;
249
250 const std::vector<PhasePresence>
251 phaseCondition() const {return this->phaseCondition_;}
252
255 void computeWaterShearVelocityFaces(const V& transi, const std::vector<ADB>& kr,
256 const std::vector<ADB>& phasePressure, const SolutionState& state,
257 std::vector<double>& water_vel, std::vector<double>& visc_mult);
258
261 void computeWaterShearVelocityWells(const SolutionState& state, WellState& xw, const ADB& cq_sw,
262 std::vector<double>& water_vel_wells, std::vector<double>& visc_mult_wells);
263
264 };
265
266
267
270 struct BlackoilPolymerSolutionState : public DefaultBlackoilSolutionState
271 {
272 explicit BlackoilPolymerSolutionState(const int np)
273 : DefaultBlackoilSolutionState(np),
274 concentration( ADB::null())
275 {
276 }
278 };
279
280
281
283 template <class Grid>
284 struct ModelTraits< BlackoilPolymerModel<Grid> >
285 {
288 typedef BlackoilModelParameters ModelParameters;
290 };
291
292} // namespace Opm
293
295
296
297#endif // OPM_BLACKOILPOLYMERMODEL_HEADER_INCLUDED
Definition: BlackoilPolymerModel.hpp:44
void computeAccum(const SolutionState &state, const int aix)
Definition: BlackoilPolymerModel_impl.hpp:225
std::vector< double > wells_bore_diameter_
Definition: BlackoilPolymerModel.hpp:149
ADB computeMc(const SolutionState &state) const
Definition: BlackoilPolymerModel_impl.hpp:551
std::vector< double > wells_rep_radius_
Definition: BlackoilPolymerModel.hpp:146
V cmax_
Definition: BlackoilPolymerModel.hpp:142
@ Concentration
Definition: BlackoilPolymerModel.hpp:133
const std::vector< PhasePresence > phaseCondition() const
Definition: BlackoilPolymerModel.hpp:251
std::vector< V > variableStateInitials(const ReservoirState &x, const WellState &xw) const
Definition: BlackoilPolymerModel_impl.hpp:163
void addWellContributionToMassBalanceEq(const std::vector< ADB > &cq_s, const SolutionState &state, WellState &xw)
Definition: BlackoilPolymerModel_impl.hpp:346
Base::ReservoirState ReservoirState
Definition: BlackoilPolymerModel.hpp:50
Base::DataBlock DataBlock
Definition: BlackoilPolymerModel.hpp:132
std::vector< double > wells_perf_length_
Definition: BlackoilPolymerModel.hpp:147
std::vector< int > variableStateIndices() const
Definition: BlackoilPolymerModel_impl.hpp:188
const bool has_shrate_
Definition: BlackoilPolymerModel.hpp:140
void prepareStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:123
BlackoilPolymerModel(const typename Base::ModelParameters &param, const Grid &grid, const BlackoilPropsAdInterface &fluid, const DerivedGeology &geo, const RockCompressibility *rock_comp_props, const PolymerPropsAd &polymer_props_ad, const Wells *wells, const NewtonIterationBlackoilInterface &linsolver, EclipseStateConstPtr eclipse_state, const bool has_disgas, const bool has_vapoil, const bool has_polymer, const bool has_plyshlog, const bool has_shrate, const std::vector< double > &wells_rep_radius, const std::vector< double > &wells_perf_length, const std::vector< double > &wells_bore_diameter, const bool terminal_output)
Definition: BlackoilPolymerModel_impl.hpp:77
Base::SolutionState SolutionState
Definition: BlackoilPolymerModel.hpp:131
SolutionState variableStateExtractVars(const ReservoirState &x, const std::vector< int > &indices, std::vector< ADB > &vars) const
Definition: BlackoilPolymerModel_impl.hpp:208
Base::WellState WellState
Definition: BlackoilPolymerModel.hpp:51
const bool has_polymer_
Definition: BlackoilPolymerModel.hpp:138
const int poly_pos_
Definition: BlackoilPolymerModel.hpp:141
void computeWaterShearVelocityWells(const SolutionState &state, WellState &xw, const ADB &cq_sw, std::vector< double > &water_vel_wells, std::vector< double > &visc_mult_wells)
Definition: BlackoilPolymerModel_impl.hpp:678
std::vector< double > shear_mult_wells_
Definition: BlackoilPolymerModel.hpp:154
const bool has_plyshlog_
Definition: BlackoilPolymerModel.hpp:139
const PolymerPropsAd & polymer_props_ad_
Definition: BlackoilPolymerModel.hpp:137
void updateState(const V &dx, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:376
BlackoilModelBase< Grid, BlackoilPolymerModel< Grid > > Base
Definition: BlackoilPolymerModel.hpp:49
void makeConstantState(SolutionState &state) const
Definition: BlackoilPolymerModel_impl.hpp:151
void assemble(const ReservoirState &reservoir_state, WellState &well_state, const bool initial_assembly)
Definition: BlackoilPolymerModel_impl.hpp:460
void computeCmax(ReservoirState &state)
Definition: BlackoilPolymerModel_impl.hpp:256
void assembleMassBalanceEq(const SolutionState &state)
Definition: BlackoilPolymerModel_impl.hpp:273
void computeWaterShearVelocityFaces(const V &transi, const std::vector< ADB > &kr, const std::vector< ADB > &phasePressure, const SolutionState &state, std::vector< double > &water_vel, std::vector< double > &visc_mult)
Definition: BlackoilPolymerModel_impl.hpp:561
std::vector< double > shear_mult_faces_
Definition: BlackoilPolymerModel.hpp:152
void afterStep(const double dt, ReservoirState &reservoir_state, WellState &well_state)
Definition: BlackoilPolymerModel_impl.hpp:138
void computeMassFlux(const int actph, const V &transi, const ADB &kr, const ADB &p, const SolutionState &state)
Definition: BlackoilPolymerModel_impl.hpp:413
Definition: PolymerBlackoilState.hpp:34
Definition: PolymerPropsAd.hpp:33
Definition: WellStateFullyImplicitBlackoilPolymer.hpp:29
Definition: CompressibleTpfaPolymer.hpp:33
Definition: BlackoilPolymerModel.hpp:271
ADB concentration
Definition: BlackoilPolymerModel.hpp:277
BlackoilPolymerSolutionState(const int np)
Definition: BlackoilPolymerModel.hpp:272
BlackoilModelParameters ModelParameters
Definition: BlackoilPolymerModel.hpp:288
BlackoilPolymerSolutionState SolutionState
Definition: BlackoilPolymerModel.hpp:289
WellStateFullyImplicitBlackoilPolymer WellState
Definition: BlackoilPolymerModel.hpp:287
PolymerBlackoilState ReservoirState
Definition: BlackoilPolymerModel.hpp:286