UpscalerBase_impl.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: UpscalerBase_impl.hpp
4//
5// Created: Thu Apr 29 10:22:06 2010
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8//
9// $Date$
10//
11// $Revision$
12//
13//===========================================================================
14
15/*
16 Copyright 2010 SINTEF ICT, Applied Mathematics.
17 Copyright 2010 Statoil ASA.
18
19 This file is part of The Open Reservoir Simulator Project (OpenRS).
20
21 OpenRS is free software: you can redistribute it and/or modify
22 it under the terms of the GNU General Public License as published by
23 the Free Software Foundation, either version 3 of the License, or
24 (at your option) any later version.
25
26 OpenRS is distributed in the hope that it will be useful,
27 but WITHOUT ANY WARRANTY; without even the implied warranty of
28 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
29 GNU General Public License for more details.
30
31 You should have received a copy of the GNU General Public License
32 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
33*/
34
35#ifndef OPM_UPSCALERBASE_IMPL_HEADER
36#define OPM_UPSCALERBASE_IMPL_HEADER
37
41
42#include <iostream>
43
44namespace Opm
45{
46
47 template <class Traits>
49 : bctype_(Fixed),
50 twodim_hack_(false),
51 residual_tolerance_(1e-8),
52 linsolver_maxit_(0),
53 linsolver_prolongate_factor_(1.0),
54 linsolver_verbosity_(0),
55 linsolver_type_(3),
56 linsolver_smooth_steps_(1),
57 gravity_(0.0)
58#if HAVE_MPI
59 , grid_(MPI_COMM_SELF)
60#endif
61 {
62 }
63
64
65
66
67 template <class Traits>
68 inline void UpscalerBase<Traits>::init(const Opm::ParameterGroup& param)
69 {
70 initImpl(param);
71 initFinal(param);
72 }
74
75
76
77 template <class Traits>
78 inline void UpscalerBase<Traits>::initImpl(const Opm::ParameterGroup& param)
79 {
80 // Request the boundary condition type parameter early since,
81 // depending on the actual type, we may have to manufacture
82 // and insert other parameters into the ParameterGroup prior
83 // to constructing the grid and associated properties.
84 //
85 int bct = param.get<int>("boundary_condition_type");
86 bctype_ = static_cast<BoundaryConditionType>(bct);
87
88 // Import control parameters pertaining to reduced physical
89 // dimensionality ("2d_hack = true" precludes periodic
90 // boundary conditions in the Z direction), and linear solves.
91 //
92 twodim_hack_ = param.getDefault("2d_hack", twodim_hack_);
93 residual_tolerance_ = param.getDefault("residual_tolerance", residual_tolerance_);
94 linsolver_verbosity_ = param.getDefault("linsolver_verbosity", linsolver_verbosity_);
95 linsolver_type_ = param.getDefault("linsolver_type", linsolver_type_);
96 linsolver_maxit_ = param.getDefault("linsolver_max_iterations", linsolver_maxit_);
97 linsolver_prolongate_factor_ = param.getDefault("linsolver_prolongate_factor", linsolver_prolongate_factor_);
98 linsolver_smooth_steps_ = param.getDefault("linsolver_smooth_steps", linsolver_smooth_steps_);
99 gravity_ = param.getDefault("gravity", gravity_);
100
101 // Ensure sufficient grid support for requested boundary
102 // condition type.
103 //
104 Opm::ParameterGroup temp_param = param;
105 if (bctype_ == Linear || bctype_ == Periodic) {
106 if (!temp_param.has("use_unique_boundary_ids")) {
107 temp_param.insertParameter("use_unique_boundary_ids", "true");
108 }
109 if (!temp_param.has("clip_z")) {
110 temp_param.insertParameter("clip_z", "true");
112 }
113 if (bctype_ == Periodic) {
114 if (!temp_param.has("periodic_extension")) {
115 temp_param.insertParameter("periodic_extension", "true");
116 }
117 }
118
119 setupGridAndProps(temp_param, grid_, res_prop_);
120 ginterf_.init(grid_);
121 }
122
124
125
126 template <class Traits>
127 inline void UpscalerBase<Traits>::initFinal(const Opm::ParameterGroup& param)
128 {
129 // Report any unused parameters.
130 std::cout << "==================== Unused parameters: ====================\n";
131 param.displayUsage();
132 std::cout << "================================================================\n";
133 }
134
135
136
137
138 template <class Traits>
139 inline void UpscalerBase<Traits>::init(const Opm::Deck& deck,
141 double perm_threshold,
142 double residual_tolerance,
143 int linsolver_verbosity,
144 int linsolver_type,
145 bool twodim_hack,
146 int linsolver_maxit,
147 double linsolver_prolongate_factor,
148 int linsolver_smooth_steps,
149 const double gravity)
150 {
151 bctype_ = bctype;
152 residual_tolerance_ = residual_tolerance;
153 linsolver_verbosity_ = linsolver_verbosity;
154 linsolver_type_ = linsolver_type;
155 linsolver_maxit_ = linsolver_maxit;
156 linsolver_prolongate_factor_ = linsolver_prolongate_factor;
157 linsolver_smooth_steps_ = linsolver_smooth_steps;
158 twodim_hack_ = twodim_hack;
159 gravity_ = gravity;
160
161
162 // Faking some parameters depending on bc type.
163 bool periodic_ext = (bctype_ == Periodic);
164 bool turn_normals = false;
165 bool clip_z = (bctype_ == Linear || bctype_ == Periodic);
166 bool unique_bids = (bctype_ == Linear || bctype_ == Periodic);
167 std::string rock_list("no_list");
169 periodic_ext, turn_normals, clip_z, unique_bids,
170 perm_threshold, rock_list,
171 useJ<ResProp>(), 1.0, 0.0,
172 grid_, res_prop_);
173 ginterf_.init(grid_);
174 }
175
176
177
178
179 template <class Traits>
180 inline const typename UpscalerBase<Traits>::GridType&
182 {
183 return grid_;
184 }
185
186
187
188
189 template <class Traits>
190 inline void
192 {
193 if ((type == Periodic && bctype_ != Periodic)
194 || (type != Periodic && bctype_ == Periodic)) {
195 OPM_THROW(std::runtime_error, "Cannot switch to or from Periodic boundary condition, "
196 "periodic must be set in init() params.");
197 } else {
198 bctype_ = type;
199 if (type == Periodic || type == Linear) {
200 grid_.setUniqueBoundaryIds(true);
201 } else {
202 grid_.setUniqueBoundaryIds(false);
203 }
204 }
205 }
206
207
208
209
210 template <class Traits>
211 inline void
213 {
214 res_prop_.permeabilityModifiable(cell_index) = k;
215 }
216
217
218
219
220 template <class Traits>
223 {
225 return upscaleEffectivePerm(fluid);
226 }
227
228
229
230
231 template <class Traits>
232 template <class FluidInterface>
235 {
236 int num_cells = ginterf_.numberOfCells();
237 // No source or sink.
238 std::vector<double> src(num_cells, 0.0);
239 // Just water.
240 std::vector<double> sat(num_cells, 1.0);
241 // Gravity.
242 Dune::FieldVector<double, 3> gravity(0.0);
243 gravity[2] = gravity_;
244
245 permtensor_t upscaled_K(3, 3, (double*)0);
246 for (int pdd = 0; pdd < Dimension; ++pdd) {
247 setupUpscalingConditions(ginterf_, bctype_, pdd, 1.0, 1.0, twodim_hack_, bcond_);
248 if (pdd == 0) {
249 // Only on first iteration, since we do not change the
250 // structure of the system, the way the flow solver is
251 // implemented.
252 flow_solver_.init(ginterf_, res_prop_, gravity, bcond_);
253 }
254
255 // Run pressure solver.
256 bool same_matrix = (bctype_ != Fixed) && (pdd != 0);
257 flow_solver_.solve(fluid, sat, bcond_, src, residual_tolerance_,
258 linsolver_verbosity_,
259 linsolver_type_, same_matrix,
260 linsolver_maxit_, linsolver_prolongate_factor_,
261 linsolver_smooth_steps_);
262 double max_mod = flow_solver_.postProcessFluxes();
263 std::cout << "Max mod = " << max_mod << std::endl;
264
265 // Compute upscaled K.
266 double Q[Dimension] = { 0 };
267 switch (bctype_) {
268 case Fixed:
269 Q[pdd] = computeAverageVelocity(flow_solver_.getSolution(), pdd, pdd);
270 break;
271 case Linear:
272 case Periodic:
273 for (int i = 0; i < Dimension; ++i) {
274 Q[i] = computeAverageVelocity(flow_solver_.getSolution(), i, pdd);
275 }
276 break;
277 default:
278 OPM_THROW(std::runtime_error, "Unknown boundary type: " + std::to_string(bctype_));
279 }
280 double delta = computeDelta(pdd);
281 for (int i = 0; i < Dimension; ++i) {
282 upscaled_K(i, pdd) = Q[i] * delta;
283 }
284 }
285 return upscaled_K;
286 }
287
288
289
290
291 template <class Traits>
292 template <class FlowSol>
293 double UpscalerBase<Traits>::computeAverageVelocity(const FlowSol& flow_solution,
294 const int flow_dir,
295 const int pdrop_dir) const
296 {
297 double side1_flux = 0.0;
298 double side2_flux = 0.0;
299 double side1_area = 0.0;
300 double side2_area = 0.0;
301
302 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
303 for (FaceIter f = c->facebegin(); f != c->faceend(); ++f) {
304 if (f->boundary()) {
305 int canon_bid = bcond_.getCanonicalBoundaryId(f->boundaryId());
306 if ((canon_bid - 1)/2 == flow_dir) {
307 double flux = flow_solution.outflux(f);
308 double area = f->area();
309 double norm_comp = f->normal()[flow_dir];
310 // std::cout << "bid " << f->boundaryId() << " area " << area << " n " << norm_comp << std::endl;
311 if (canon_bid - 1 == 2*flow_dir) {
312 if (flow_dir == pdrop_dir && flux > 0.0) {
313#ifdef VERBOSE
314 std::cerr << "Flow may be in wrong direction at bid: " << f->boundaryId()<<" (canonical: "<<canon_bid
315 << ") Magnitude: " << std::fabs(flux) << std::endl;
316#endif
317 // OPM_THROW(std::runtime_error, "Detected outflow at entry face: " << face);
318 }
319 side1_flux += flux*norm_comp;
320 side1_area += area;
321 } else {
322 assert(canon_bid - 1 == 2*flow_dir + 1);
323 if (flow_dir == pdrop_dir && flux < 0.0) {
324#ifdef VERBOSE
325 std::cerr << "Flow may be in wrong direction at bid: " << f->boundaryId()
326 << " Magnitude: " << std::fabs(flux) << std::endl;
327#endif
328 // OPM_THROW(std::runtime_error, "Detected inflow at exit face: " << face);
329 }
330 side2_flux += flux*norm_comp;
331 side2_area += area;
332 }
333 }
334 }
335 }
336 }
337 // q is the average velocity.
338 return 0.5*(side1_flux/side1_area + side2_flux/side2_area);
339 }
340
341
342
343
344 template <class Traits>
345 inline double UpscalerBase<Traits>::computeDelta(const int flow_dir) const
346 {
347 double side1_pos = 0.0;
348 double side2_pos = 0.0;
349 double side1_area = 0.0;
350 double side2_area = 0.0;
351 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
352 for (FaceIter f = c->facebegin(); f != c->faceend(); ++f) {
353 if (f->boundary()) {
354 int canon_bid = bcond_.getCanonicalBoundaryId(f->boundaryId());
355 if ((canon_bid - 1)/2 == flow_dir) {
356 double area = f->area();
357 double pos_comp = f->centroid()[flow_dir];
358 if (canon_bid - 1 == 2*flow_dir) {
359 side1_pos += area*pos_comp;
360 side1_area += area;
361 } else {
362 side2_pos += area*pos_comp;
363 side2_area += area;
364 }
365 }
366 }
367 }
368 }
369 // delta is the average length.
370 return side2_pos/side2_area - side1_pos/side1_area;
371 }
372
373
374
375
376 template <class Traits>
378 {
379 double total_vol = 0.0;
380 double total_pore_vol = 0.0;
381 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
382 total_vol += c->volume();
383 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
384 }
385 return total_pore_vol/total_vol;
386 }
387
388
389 template <class Traits>
391 {
392 double total_net_vol = 0.0;
393 double total_pore_vol = 0.0;
394 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
395 total_net_vol += c->volume()*res_prop_.ntg(c->index());
396 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
397 }
398 if (total_net_vol>0.0) return total_pore_vol/total_net_vol;
399 else return 0.0;
400 }
401
402 template <class Traits>
404 {
405 double total_vol = 0.0;
406 double total_net_vol = 0.0;
407 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
408 total_vol += c->volume();
409 total_net_vol += c->volume()*res_prop_.ntg(c->index());
410 }
411 return total_net_vol/total_vol;
412 }
413
414 template <class Traits>
415 double UpscalerBase<Traits>::upscaleSWCR(const bool NTG) const
416 {
417 double total_swcr = 0.0;
418 double total_pore_vol = 0.0;
419 if (NTG) {
420 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
421 total_swcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index())*res_prop_.swcr(c->index());
422 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
423 }
424 }
425 else {
426 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
427 total_swcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.swcr(c->index());
428 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
429 }
430 }
431 return total_swcr/total_pore_vol;
432 }
433
434 template <class Traits>
435 double UpscalerBase<Traits>::upscaleSOWCR(const bool NTG) const
436 {
437 double total_sowcr = 0.0;
438 double total_pore_vol = 0.0;
439 if (NTG) {
440 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
441 total_sowcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index())*res_prop_.sowcr(c->index());
442 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
443 }
444 }
445 else {
446 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
447 total_sowcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.sowcr(c->index());
448 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
449 }
450 }
451 return total_sowcr/total_pore_vol;
452 }
453
454} // namespace Opm
455
456
457
458#endif // OPM_UPSCALERBASE_IMPL_HEADER
Definition: GridInterfaceEuler.hpp:368
Intersection (face) iterator for solver-near grid interface.
Definition: GridInterfaceEuler.hpp:241
Definition: ReservoirPropertyTracerFluid.hpp:40
double computeDelta(const int flow_dir) const
Definition: UpscalerBase_impl.hpp:345
permtensor_t upscaleSinglePhase()
Definition: UpscalerBase_impl.hpp:222
void setPermeability(const int cell_index, const permtensor_t &k)
Definition: UpscalerBase_impl.hpp:212
double computeAverageVelocity(const FlowSol &flow_solution, const int flow_dir, const int pdrop_dir) const
Definition: UpscalerBase_impl.hpp:293
double upscalePorosity() const
Definition: UpscalerBase_impl.hpp:377
ResProp::MutablePermTensor permtensor_t
A type for the upscaled permeability.
Definition: UpscalerBase.hpp:66
const GridType & grid() const
Access the grid.
Definition: UpscalerBase_impl.hpp:181
double upscaleNTG() const
Definition: UpscalerBase_impl.hpp:403
virtual void initFinal(const Opm::ParameterGroup &param)
Definition: UpscalerBase_impl.hpp:127
BoundaryConditionType
Definition: UpscalerBase.hpp:68
permtensor_t upscaleEffectivePerm(const FluidInterface &fluid)
Definition: UpscalerBase_impl.hpp:234
Dune::CpGrid GridType
Definition: UpscalerBase.hpp:60
double upscaleSOWCR(const bool NTG) const
Definition: UpscalerBase_impl.hpp:435
virtual void initImpl(const Opm::ParameterGroup &param)
Definition: UpscalerBase_impl.hpp:78
void init(const Opm::ParameterGroup &param)
Initializes the upscaler from parameters.
Definition: UpscalerBase_impl.hpp:68
void setBoundaryConditionType(BoundaryConditionType type)
Definition: UpscalerBase_impl.hpp:191
double upscaleNetPorosity() const
Definition: UpscalerBase_impl.hpp:390
double upscaleSWCR(const bool NTG) const
Definition: UpscalerBase_impl.hpp:415
UpscalerBase()
Default constructor.
Definition: UpscalerBase_impl.hpp:48
if(loadcase > -1)
Definition: elasticity_upscale_impl.hpp:448
auto deck
Definition: elasticity_upscale_impl.hpp:590
Definition: ImplicitAssembly.hpp:43
void setupGridAndPropsEclipse(const Opm::Deck &deck, bool periodic_extension, bool turn_normals, bool clip_z, bool unique_bids, double perm_threshold, const std::string &rock_list, bool use_jfunction_scaling, double sigma, double theta, Dune::CpGrid &grid, ResProp< 3 > &res_prop)
Definition: setupGridAndProps.hpp:152
void setupUpscalingConditions(const GridInterface &g, int bct, int pddir, double pdrop, double bdy_sat, bool twodim_hack, BCs &bcs)
Definition: setupBoundaryConditions.hpp:99
void setupGridAndProps(const Opm::ParameterGroup &param, Dune::CpGrid &grid, ResProp< 3 > &res_prop)
Definition: setupGridAndProps.hpp:71