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 {
59 }
60
61
62
63
64 template <class Traits>
65 inline void UpscalerBase<Traits>::init(const Opm::ParameterGroup& param)
66 {
67 initImpl(param);
68 initFinal(param);
69 }
70
71
72
74 template <class Traits>
75 inline void UpscalerBase<Traits>::initImpl(const Opm::ParameterGroup& param)
76 {
77 // Request the boundary condition type parameter early since,
78 // depending on the actual type, we may have to manufacture
79 // and insert other parameters into the ParameterGroup prior
80 // to constructing the grid and associated properties.
81 //
82 int bct = param.get<int>("boundary_condition_type");
83 bctype_ = static_cast<BoundaryConditionType>(bct);
84
85 // Import control parameters pertaining to reduced physical
86 // dimensionality ("2d_hack = true" precludes periodic
87 // boundary conditions in the Z direction), and linear solves.
88 //
89 twodim_hack_ = param.getDefault("2d_hack", twodim_hack_);
90 residual_tolerance_ = param.getDefault("residual_tolerance", residual_tolerance_);
91 linsolver_verbosity_ = param.getDefault("linsolver_verbosity", linsolver_verbosity_);
92 linsolver_type_ = param.getDefault("linsolver_type", linsolver_type_);
93 linsolver_maxit_ = param.getDefault("linsolver_max_iterations", linsolver_maxit_);
94 linsolver_prolongate_factor_ = param.getDefault("linsolver_prolongate_factor", linsolver_prolongate_factor_);
95 linsolver_smooth_steps_ = param.getDefault("linsolver_smooth_steps", linsolver_smooth_steps_);
96 gravity_ = param.getDefault("gravity", gravity_);
97
98 // Ensure sufficient grid support for requested boundary
99 // condition type.
100 //
101 Opm::ParameterGroup temp_param = param;
102 if (bctype_ == Linear || bctype_ == Periodic) {
103 if (!temp_param.has("use_unique_boundary_ids")) {
104 temp_param.insertParameter("use_unique_boundary_ids", "true");
105 }
106 if (!temp_param.has("clip_z")) {
107 temp_param.insertParameter("clip_z", "true");
108 }
109 }
110 if (bctype_ == Periodic) {
111 if (!temp_param.has("periodic_extension")) {
112 temp_param.insertParameter("periodic_extension", "true");
113 }
114 }
116 setupGridAndProps(temp_param, grid_, res_prop_);
117 ginterf_.init(grid_);
118 }
120
121
122
123 template <class Traits>
124 inline void UpscalerBase<Traits>::initFinal(const Opm::ParameterGroup& param)
125 {
126 // Report any unused parameters.
127 std::cout << "==================== Unused parameters: ====================\n";
128 param.displayUsage();
129 std::cout << "================================================================\n";
130 }
131
132
133
134
135 template <class Traits>
136 inline void UpscalerBase<Traits>::init(const Opm::Deck& deck,
138 double perm_threshold,
139 double residual_tolerance,
140 int linsolver_verbosity,
141 int linsolver_type,
142 bool twodim_hack,
143 int linsolver_maxit,
144 double linsolver_prolongate_factor,
145 int linsolver_smooth_steps,
146 const double gravity)
148 bctype_ = bctype;
149 residual_tolerance_ = residual_tolerance;
150 linsolver_verbosity_ = linsolver_verbosity;
151 linsolver_type_ = linsolver_type;
152 linsolver_maxit_ = linsolver_maxit;
153 linsolver_prolongate_factor_ = linsolver_prolongate_factor;
154 linsolver_smooth_steps_ = linsolver_smooth_steps;
155 twodim_hack_ = twodim_hack;
156 gravity_ = gravity;
157
158
159 // Faking some parameters depending on bc type.
160 bool periodic_ext = (bctype_ == Periodic);
161 bool turn_normals = false;
162 bool clip_z = (bctype_ == Linear || bctype_ == Periodic);
163 bool unique_bids = (bctype_ == Linear || bctype_ == Periodic);
164 std::string rock_list("no_list");
166 periodic_ext, turn_normals, clip_z, unique_bids,
167 perm_threshold, rock_list,
168 useJ<ResProp>(), 1.0, 0.0,
169 grid_, res_prop_);
170 ginterf_.init(grid_);
171 }
172
173
174
175
176 template <class Traits>
177 inline const typename UpscalerBase<Traits>::GridType&
179 {
180 return grid_;
181 }
182
183
184
185
186 template <class Traits>
187 inline void
189 {
190 if ((type == Periodic && bctype_ != Periodic)
191 || (type != Periodic && bctype_ == Periodic)) {
192 OPM_THROW(std::runtime_error, "Cannot switch to or from Periodic boundary condition, "
193 "periodic must be set in init() params.");
194 } else {
195 bctype_ = type;
196 if (type == Periodic || type == Linear) {
197 grid_.setUniqueBoundaryIds(true);
198 } else {
199 grid_.setUniqueBoundaryIds(false);
200 }
201 }
202 }
203
204
205
206
207 template <class Traits>
208 inline void
210 {
211 res_prop_.permeabilityModifiable(cell_index) = k;
212 }
213
214
215
216
217 template <class Traits>
220 {
222 return upscaleEffectivePerm(fluid);
223 }
224
225
226
227
228 template <class Traits>
229 template <class FluidInterface>
232 {
233 int num_cells = ginterf_.numberOfCells();
234 // No source or sink.
235 std::vector<double> src(num_cells, 0.0);
236 // Just water.
237 std::vector<double> sat(num_cells, 1.0);
238 // Gravity.
239 Dune::FieldVector<double, 3> gravity(0.0);
240 gravity[2] = gravity_;
241
242 permtensor_t upscaled_K(3, 3, (double*)0);
243 for (int pdd = 0; pdd < Dimension; ++pdd) {
244 setupUpscalingConditions(ginterf_, bctype_, pdd, 1.0, 1.0, twodim_hack_, bcond_);
245 if (pdd == 0) {
246 // Only on first iteration, since we do not change the
247 // structure of the system, the way the flow solver is
248 // implemented.
249 flow_solver_.init(ginterf_, res_prop_, gravity, bcond_);
250 }
251
252 // Run pressure solver.
253 bool same_matrix = (bctype_ != Fixed) && (pdd != 0);
254 flow_solver_.solve(fluid, sat, bcond_, src, residual_tolerance_,
255 linsolver_verbosity_,
256 linsolver_type_, same_matrix,
257 linsolver_maxit_, linsolver_prolongate_factor_,
258 linsolver_smooth_steps_);
259 double max_mod = flow_solver_.postProcessFluxes();
260 std::cout << "Max mod = " << max_mod << std::endl;
261
262 // Compute upscaled K.
263 double Q[Dimension] = { 0 };
264 switch (bctype_) {
265 case Fixed:
266 Q[pdd] = computeAverageVelocity(flow_solver_.getSolution(), pdd, pdd);
267 break;
268 case Linear:
269 case Periodic:
270 for (int i = 0; i < Dimension; ++i) {
271 Q[i] = computeAverageVelocity(flow_solver_.getSolution(), i, pdd);
272 }
273 break;
274 default:
275 OPM_THROW(std::runtime_error, "Unknown boundary type: " + std::to_string(bctype_));
276 }
277 double delta = computeDelta(pdd);
278 for (int i = 0; i < Dimension; ++i) {
279 upscaled_K(i, pdd) = Q[i] * delta;
280 }
281 }
282 return upscaled_K;
283 }
284
285
286
287
288 template <class Traits>
289 template <class FlowSol>
290 double UpscalerBase<Traits>::computeAverageVelocity(const FlowSol& flow_solution,
291 const int flow_dir,
292 const int pdrop_dir) const
293 {
294 double side1_flux = 0.0;
295 double side2_flux = 0.0;
296 double side1_area = 0.0;
297 double side2_area = 0.0;
298
299 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
300 for (FaceIter f = c->facebegin(); f != c->faceend(); ++f) {
301 if (f->boundary()) {
302 int canon_bid = bcond_.getCanonicalBoundaryId(f->boundaryId());
303 if ((canon_bid - 1)/2 == flow_dir) {
304 double flux = flow_solution.outflux(f);
305 double area = f->area();
306 double norm_comp = f->normal()[flow_dir];
307 // std::cout << "bid " << f->boundaryId() << " area " << area << " n " << norm_comp << std::endl;
308 if (canon_bid - 1 == 2*flow_dir) {
309 if (flow_dir == pdrop_dir && flux > 0.0) {
310#ifdef VERBOSE
311 std::cerr << "Flow may be in wrong direction at bid: " << f->boundaryId()<<" (canonical: "<<canon_bid
312 << ") Magnitude: " << std::fabs(flux) << std::endl;
313#endif
314 // OPM_THROW(std::runtime_error, "Detected outflow at entry face: " << face);
315 }
316 side1_flux += flux*norm_comp;
317 side1_area += area;
318 } else {
319 assert(canon_bid - 1 == 2*flow_dir + 1);
320 if (flow_dir == pdrop_dir && flux < 0.0) {
321#ifdef VERBOSE
322 std::cerr << "Flow may be in wrong direction at bid: " << f->boundaryId()
323 << " Magnitude: " << std::fabs(flux) << std::endl;
324#endif
325 // OPM_THROW(std::runtime_error, "Detected inflow at exit face: " << face);
326 }
327 side2_flux += flux*norm_comp;
328 side2_area += area;
329 }
330 }
331 }
332 }
333 }
334 // q is the average velocity.
335 return 0.5*(side1_flux/side1_area + side2_flux/side2_area);
336 }
337
338
339
340
341 template <class Traits>
342 inline double UpscalerBase<Traits>::computeDelta(const int flow_dir) const
343 {
344 double side1_pos = 0.0;
345 double side2_pos = 0.0;
346 double side1_area = 0.0;
347 double side2_area = 0.0;
348 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
349 for (FaceIter f = c->facebegin(); f != c->faceend(); ++f) {
350 if (f->boundary()) {
351 int canon_bid = bcond_.getCanonicalBoundaryId(f->boundaryId());
352 if ((canon_bid - 1)/2 == flow_dir) {
353 double area = f->area();
354 double pos_comp = f->centroid()[flow_dir];
355 if (canon_bid - 1 == 2*flow_dir) {
356 side1_pos += area*pos_comp;
357 side1_area += area;
358 } else {
359 side2_pos += area*pos_comp;
360 side2_area += area;
361 }
362 }
363 }
364 }
365 }
366 // delta is the average length.
367 return side2_pos/side2_area - side1_pos/side1_area;
368 }
369
370
371
372
373 template <class Traits>
375 {
376 double total_vol = 0.0;
377 double total_pore_vol = 0.0;
378 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
379 total_vol += c->volume();
380 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
381 }
382 return total_pore_vol/total_vol;
383 }
384
385
386 template <class Traits>
388 {
389 double total_net_vol = 0.0;
390 double total_pore_vol = 0.0;
391 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
392 total_net_vol += c->volume()*res_prop_.ntg(c->index());
393 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
394 }
395 if (total_net_vol>0.0) return total_pore_vol/total_net_vol;
396 else return 0.0;
397 }
398
399 template <class Traits>
401 {
402 double total_vol = 0.0;
403 double total_net_vol = 0.0;
404 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
405 total_vol += c->volume();
406 total_net_vol += c->volume()*res_prop_.ntg(c->index());
407 }
408 return total_net_vol/total_vol;
409 }
410
411 template <class Traits>
412 double UpscalerBase<Traits>::upscaleSWCR(const bool NTG) const
413 {
414 double total_swcr = 0.0;
415 double total_pore_vol = 0.0;
416 if (NTG) {
417 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
418 total_swcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index())*res_prop_.swcr(c->index());
419 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
420 }
421 }
422 else {
423 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
424 total_swcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.swcr(c->index());
425 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
426 }
427 }
428 return total_swcr/total_pore_vol;
429 }
430
431 template <class Traits>
432 double UpscalerBase<Traits>::upscaleSOWCR(const bool NTG) const
433 {
434 double total_sowcr = 0.0;
435 double total_pore_vol = 0.0;
436 if (NTG) {
437 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
438 total_sowcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index())*res_prop_.sowcr(c->index());
439 total_pore_vol += c->volume()*res_prop_.porosity(c->index())*res_prop_.ntg(c->index());
440 }
441 }
442 else {
443 for (CellIter c = ginterf_.cellbegin(); c != ginterf_.cellend(); ++c) {
444 total_sowcr += c->volume()*res_prop_.porosity(c->index())*res_prop_.sowcr(c->index());
445 total_pore_vol += c->volume()*res_prop_.porosity(c->index());
446 }
447 }
448 return total_sowcr/total_pore_vol;
449 }
450
451} // namespace Opm
452
453
454
455#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:342
permtensor_t upscaleSinglePhase()
Definition: UpscalerBase_impl.hpp:219
void setPermeability(const int cell_index, const permtensor_t &k)
Definition: UpscalerBase_impl.hpp:209
double computeAverageVelocity(const FlowSol &flow_solution, const int flow_dir, const int pdrop_dir) const
Definition: UpscalerBase_impl.hpp:290
double upscalePorosity() const
Definition: UpscalerBase_impl.hpp:374
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:178
double upscaleNTG() const
Definition: UpscalerBase_impl.hpp:400
virtual void initFinal(const Opm::ParameterGroup &param)
Definition: UpscalerBase_impl.hpp:124
BoundaryConditionType
Definition: UpscalerBase.hpp:68
permtensor_t upscaleEffectivePerm(const FluidInterface &fluid)
Definition: UpscalerBase_impl.hpp:231
Dune::CpGrid GridType
Definition: UpscalerBase.hpp:60
double upscaleSOWCR(const bool NTG) const
Definition: UpscalerBase_impl.hpp:432
virtual void initImpl(const Opm::ParameterGroup &param)
Definition: UpscalerBase_impl.hpp:75
void init(const Opm::ParameterGroup &param)
Initializes the upscaler from parameters.
Definition: UpscalerBase_impl.hpp:65
void setBoundaryConditionType(BoundaryConditionType type)
Definition: UpscalerBase_impl.hpp:188
double upscaleNetPorosity() const
Definition: UpscalerBase_impl.hpp:387
double upscaleSWCR(const bool NTG) const
Definition: UpscalerBase_impl.hpp:412
UpscalerBase()
Default constructor.
Definition: UpscalerBase_impl.hpp:48
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