EulerUpstream_impl.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: EulerUpstream_impl.hpp
4//
5// Created: Tue Jun 16 14:25:24 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// B�rd Skaflestad <bard.skaflestad@sintef.no>
9// Halvor M Nilsen <hnil@sintef.no>
10//
11// $Date$
12//
13// $Revision$
14//
15//===========================================================================
16
17/*
18 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
19 Copyright 2009, 2010 Statoil ASA.
20
21 This file is part of The Open Reservoir Simulator Project (OpenRS).
22
23 OpenRS is free software: you can redistribute it and/or modify
24 it under the terms of the GNU General Public License as published by
25 the Free Software Foundation, either version 3 of the License, or
26 (at your option) any later version.
27
28 OpenRS is distributed in the hope that it will be useful,
29 but WITHOUT ANY WARRANTY; without even the implied warranty of
30 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
31 GNU General Public License for more details.
32
33 You should have received a copy of the GNU General Public License
34 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
35*/
36
37#ifndef OPENRS_EULERUPSTREAM_IMPL_HEADER
38#define OPENRS_EULERUPSTREAM_IMPL_HEADER
39
40
41
42#include <opm/common/ErrorMacros.hpp>
43#include <opm/core/utility/Average.hpp>
44#include <opm/core/utility/Units.hpp>
45#include <dune/grid/common/Volumes.hpp>
47#include <opm/core/utility/StopWatch.hpp>
48
49#include <cassert>
50#include <cmath>
51#include <algorithm>
52#include <limits>
53#include <iostream>
54
55namespace Opm
56{
57
58
59 template <class GI, class RP, class BC>
61 : method_viscous_(true),
62 method_gravity_(true),
63 method_capillary_(true),
64 use_cfl_viscous_(true),
65 use_cfl_gravity_(true),
66 use_cfl_capillary_(true),
67 courant_number_(0.5),
68 minimum_small_steps_(1),
69 maximum_small_steps_(10000),
70 check_sat_(true),
71 clamp_sat_(false)
72 {
73 }
74
75 template <class GI, class RP, class BC>
76 inline EulerUpstream<GI, RP, BC>::EulerUpstream(const GI& g, const RP& r, const BC& b)
77 : method_viscous_(true),
78 method_gravity_(true),
79 method_capillary_(true),
80 use_cfl_viscous_(true),
81 use_cfl_gravity_(true),
82 use_cfl_capillary_(true),
83 courant_number_(0.5),
84 minimum_small_steps_(1),
85 maximum_small_steps_(10000),
86 check_sat_(true),
87 clamp_sat_(false)
88 {
89 residual_computer_.initObj(g, r, b);
90 }
91
92
93
94 template <class GI, class RP, class BC>
95 inline void EulerUpstream<GI, RP, BC>::init(const Opm::parameter::ParameterGroup& param)
96 {
97 courant_number_ = param.getDefault("courant_number", courant_number_);
98 method_viscous_ = param.getDefault("method_viscous", method_viscous_);
99 method_gravity_ = param.getDefault("method_gravity", method_gravity_);
100 method_capillary_ = param.getDefault("method_capillary", method_capillary_);
101 use_cfl_viscous_ = param.getDefault("use_cfl_viscous", use_cfl_viscous_);
102 use_cfl_gravity_ = param.getDefault("use_cfl_gravity", use_cfl_gravity_);
103 use_cfl_capillary_ = param.getDefault("use_cfl_capillary", use_cfl_capillary_);
104 minimum_small_steps_ = param.getDefault("minimum_small_steps", minimum_small_steps_);
105 maximum_small_steps_ = param.getDefault("maximum_small_steps", maximum_small_steps_);
106 check_sat_ = param.getDefault("check_sat", check_sat_);
107 clamp_sat_ = param.getDefault("clamp_sat", clamp_sat_);
108 }
109
110 template <class GI, class RP, class BC>
111 inline void EulerUpstream<GI, RP, BC>::init(const Opm::parameter::ParameterGroup& param,
112 const GI& g, const RP& r, const BC& b)
113 {
114 init(param);
115 initObj(g, r, b);
116 }
117
118
119 template <class GI, class RP, class BC>
120 inline void EulerUpstream<GI, RP, BC>::initObj(const GI& g, const RP& r, const BC& b)
121 {
122 residual_computer_.initObj(g, r, b);
123 porevol_.resize(g.numberOfCells());
124 for (CIt c = g.cellbegin(); c != g.cellend(); ++c) {
125 porevol_[c->index()] = c->volume()*r.porosity(c->index());
126 }
127 }
128
129
130
131 template <class GI, class RP, class BC>
133 {
134 using namespace std;
135 cout << endl;
136 cout <<"Displaying some members of EulerUpstream" << endl;
137 cout << endl;
138 cout << "courant_number = " << courant_number_ << endl;
139 }
140
141
142
143 template <class GI, class RP, class BC>
145 {
146 courant_number_ = cn;
147 }
148
149
150
151 template <class GI, class RP, class BC>
152 template <class PressureSolution>
153 void EulerUpstream<GI, RP, BC>::transportSolve(std::vector<double>& saturation,
154 const double time,
155 const typename GI::Vector& gravity,
156 const PressureSolution& pressure_sol,
157 const Opm::SparseVector<double>& injection_rates) const
158 {
159 // Compute the cfl time-step.
160 double cfl_dt = computeCflTime(saturation, time, gravity, pressure_sol);
161
162 // Compute the number of small steps to take, and the actual small timestep.
163 int nr_transport_steps;
164 if (cfl_dt > time){
165 nr_transport_steps = minimum_small_steps_;
166 } else {
167 double steps = std::min<double>(std::ceil(time/cfl_dt), std::numeric_limits<int>::max());
168 nr_transport_steps = int(steps);
169 nr_transport_steps = std::max(nr_transport_steps, minimum_small_steps_);
170 nr_transport_steps = std::min(nr_transport_steps, maximum_small_steps_);
171 }
172 double dt_transport = time/nr_transport_steps;
173
174 // Do the timestepping. The try-catch blocks are there to handle
175 // the situation that smallTimeStep throws, which may happen due
176 // to saturation out of bounds (if check_sat_ is true).
177 // We cannot guarantee that this does not happen, since we do not
178 // (yet) compute a capillary cfl condition.
179 // Using exception for "alternate control flow" like this is bad
180 // design, should rather use error return values for this.
181 std::vector<double> saturation_initial(saturation);
182 bool finished = false;
183 int repeats = 0;
184 const int max_repeats = 10;
185 Opm::time::StopWatch clock;
186 clock.start();
187 while (!finished) {
188 try {
189#ifdef VERBOSE
190 std::cout << "Doing " << nr_transport_steps
191 << " steps for saturation equation with stepsize "
192 << dt_transport << " in seconds." << std::endl;
193#endif // VERBOSE
194 for (int q = 0; q < nr_transport_steps; ++q) {
195 smallTimeStep(saturation,
196 dt_transport,
197 gravity,
198 pressure_sol,
199 injection_rates);
200 }
201 finished = true;
202 }
203 catch (...) {
204 ++repeats;
205 if (repeats > max_repeats) {
206 throw;
207 }
208 OPM_MESSAGE("Warning: Transport failed, retrying with more steps.");
209 nr_transport_steps *= 2;
210 dt_transport = time/nr_transport_steps;
211 saturation = saturation_initial;
212 }
213 }
214 clock.stop();
215#ifdef VERBOSE
216 std::cout << "Seconds taken by transport solver: " << clock.secsSinceStart() << std::endl;
217#endif // VERBOSE
218 }
219
220
221
222
223
224 /*
225
226
227 template <class GI, class RP, class BC>
228 inline void EulerUpstream<GI, RP, BC>::initFinal()
229 {
230 // Build bid_to_face_ mapping for handling periodic conditions.
231 int maxbid = 0;
232 for (typename GI::CellIterator c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
233 for (typename GI::CellIterator::FaceIterator f = c->facebegin(); f != c->faceend(); ++f) {
234 int bid = f->boundaryId();
235 maxbid = std::max(maxbid, bid);
236 }
237 }
238 bid_to_face_.clear();
239 bid_to_face_.resize(maxbid + 1);
240 for (typename GI::CellIterator c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c) {
241 for (typename GI::CellIterator::FaceIterator f = c->facebegin(); f != c->faceend(); ++f) {
242 if (f->boundary() && pboundary_->satCond(*f).isPeriodic()) {
243 bid_to_face_[f->boundaryId()] = f;
244 }
245 }
246 }
247
248 // Build cell_iters_.
249 const int num_cells_per_iter = std::min(50, pgrid_->numberOfCells());
250 int counter = 0;
251 for (typename GI::CellIterator c = pgrid_->cellbegin(); c != pgrid_->cellend(); ++c, ++counter) {
252 if (counter % num_cells_per_iter == 0) {
253 cell_iters_.push_back(c);
254 }
255 }
256 cell_iters_.push_back(pgrid_->cellend());
257 }
258
259 */
260
261
262
263 template <class GI, class RP, class BC>
264 template <class PressureSolution>
265 inline double EulerUpstream<GI, RP, BC>::computeCflTime(const std::vector<double>& /*saturation*/,
266#ifdef VERBOSE
267 const double time,
268#else
269 const double,
270#endif
271 const typename GI::Vector& gravity,
272 const PressureSolution& pressure_sol) const
273 {
274 // Deal with cfl stuff, compute the necessary number of small time steps.
275 double cfl_dt_v = 1e99;
276 double cfl_dt_g = 1e99;
277 double cfl_dt_c = 1e99;
278
279 // Viscous cfl.
280 if (method_viscous_ && use_cfl_viscous_) {
281 cfl_dt_v = cfl_calculator::findCFLtimeVelocity(residual_computer_.grid(),
282 residual_computer_.reservoirProperties(),
283 pressure_sol);
284#ifdef VERBOSE
285 std::cout << "CFL dt for velocity is "
286 << cfl_dt_v << " seconds ("
287 << Opm::unit::convert::to(cfl_dt_v, Opm::unit::day)
288 << " days)." << std::endl;
289#endif // VERBOSE
290 }
291
292 // Gravity cfl.
293 if (method_gravity_ && use_cfl_gravity_) {
294 cfl_dt_g = cfl_calculator::findCFLtimeGravity(residual_computer_.grid(),
295 residual_computer_.reservoirProperties(),
296 gravity);
297#ifdef VERBOSE
298 std::cout << "CFL dt for gravity is "
299 << cfl_dt_g << " seconds ("
300 << Opm::unit::convert::to(cfl_dt_g, Opm::unit::day)
301 << " days)." << std::endl;
302#endif // VERBOSE
303 }
304
305 // Capillary cfl.
306 if (method_capillary_ && use_cfl_capillary_) {
307 cfl_dt_c = cfl_calculator::findCFLtimeCapillary(residual_computer_.grid(),
308 residual_computer_.reservoirProperties());
309#ifdef VERBOSE
310 std::cout << "CFL dt for capillary term is "
311 << cfl_dt_c << " seconds ("
312 << Opm::unit::convert::to(cfl_dt_c, Opm::unit::day)
313 << " days)." << std::endl;
314#endif // VERBOSE
315 }
316
317 double cfl_dt = std::min(std::min(cfl_dt_v, cfl_dt_g), cfl_dt_c);
318 cfl_dt *= courant_number_;
319#ifdef VERBOSE
320 std::cout << "Total impes time is "
321 << time << " seconds ("
322 << Opm::unit::convert::to(time, Opm::unit::day)
323 << " days)." << std::endl;
324
325 std::cout << "Final modified CFL dt is "
326 << cfl_dt << " seconds ("
327 << Opm::unit::convert::to(cfl_dt, Opm::unit::day)
328 << " days)." << std::endl;
329#endif // VERBOSE
330 return cfl_dt;
331 }
332
333
334
335
336 template <class GI, class RP, class BC>
337 inline void EulerUpstream<GI, RP, BC>::checkAndPossiblyClampSat(std::vector<double>& s) const
338 {
339 int num_cells = s.size();
340 for (int cell = 0; cell < num_cells; ++cell) {
341 if (s[cell] > 1.0 || s[cell] < 0.0) {
342 if (clamp_sat_) {
343 s[cell] = std::max(std::min(s[cell], 1.0), 0.0);
344 } else if (s[cell] > 1.001 || s[cell] < -0.001) {
345 OPM_THROW(std::runtime_error, "Saturation out of range in EulerUpstream: Cell " << cell << " sat " << s[cell]);
346 }
347 }
348 }
349 }
350
351
352
353
354
355 template <class GI, class RP, class BC>
356 template <class PressureSolution>
357 inline void EulerUpstream<GI, RP, BC>::smallTimeStep(std::vector<double>& saturation,
358 const double dt,
359 const typename GI::Vector& gravity,
360 const PressureSolution& pressure_sol,
361 const Opm::SparseVector<double>& injection_rates) const
362 {
363 if (method_capillary_) {
364 residual_computer_.computeCapPressures(saturation);
365 }
366 residual_computer_.computeResidual(saturation, gravity, pressure_sol, injection_rates,
367 method_viscous_, method_gravity_, method_capillary_,
368 residual_);
369// double max_ok_dt = 1e100;
370// const double tol = 1e-10;
371 int num_cells = saturation.size();
372 for (int i = 0; i < num_cells; ++i) {
373 const double sat_change = dt*residual_[i]/porevol_[i];
374 saturation[i] += sat_change;
375// if (sat_change > tol) {
376// max_ok_dt = std::min(max_ok_dt, (1.0 - saturation[i])/sc);
377// } else if (sat_change < -tol) {
378// max_ok_dt = std::min(max_ok_dt, -saturation[i]/sc);
379// }
380 }
381// std::cout << "Maximum nonviolating timestep is " << max_ok_dt << " seconds\n";
382 if (check_sat_ || clamp_sat_) {
383 checkAndPossiblyClampSat(saturation);
384 }
385 }
386
387
388} // end namespace Opm
389
390
391#endif // OPENRS_EULERUPSTREAM_IMPL_HEADER
void transportSolve(std::vector< double > &saturation, const double time, const typename GridInterface::Vector &gravity, const PressureSolution &pressure_sol, const Opm::SparseVector< double > &injection_rates) const
Solve transport equation, evolving.
void setCourantNumber(double cn)
Set the Courant number. That is dt = dt_cfl*courant_number. For this explicit method it should be < 1...
Definition: EulerUpstream_impl.hpp:144
double computeCflTime(const std::vector< double > &saturation, const double time, const typename GridInterface::Vector &gravity, const PressureSolution &pressure_sol) const
void initObj(const GridInterface &grid, const ReservoirProperties &resprop, const BoundaryConditions &boundary)
Definition: EulerUpstream_impl.hpp:120
EulerUpstream()
Definition: EulerUpstream_impl.hpp:60
void init(const Opm::parameter::ParameterGroup &param)
Definition: EulerUpstream_impl.hpp:95
GridInterface::CellIterator CIt
Definition: EulerUpstream.hpp:105
void display()
Definition: EulerUpstream_impl.hpp:132
void smallTimeStep(std::vector< double > &saturation, const double time, const typename GridInterface::Vector &gravity, const PressureSolution &pressure_sol, const Opm::SparseVector< double > &injection_rates) const
EulerUpstreamResidual< GridInterface, ReservoirProperties, BoundaryConditions > residual_computer_
Definition: EulerUpstream.hpp:128
void checkAndPossiblyClampSat(std::vector< double > &s) const
Definition: EulerUpstream_impl.hpp:337
double findCFLtimeVelocity(const Grid &grid, const ReservoirProperties &resprop, const PressureSolution &pressure_sol)
Definition: CflCalculator.hpp:55
double findCFLtimeGravity(const Grid &grid, const ReservoirProperties &resprop, const typename Grid::Vector &gravity)
Definition: CflCalculator.hpp:91
double findCFLtimeCapillary(const Grid &grid, const ReservoirProperties &resprop)
Definition: CflCalculator.hpp:143
Definition: BlackoilFluid.hpp:32