CflCalculator.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: CflCalculator.hpp
4//
5// Created: Wed Jul 1 09:35:59 2009
6//
7// Author(s): Halvor M Nilsen <hnil@sintef.no>
8// Atgeirr F Rasmussen <atgeirr@sintef.no>
9//
10// $Date$
11//
12// $Revision$
13//
14//===========================================================================
15
16/*
17 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18 Copyright 2009, 2010 Statoil ASA.
19
20 This file is part of The Open Reservoir Simulator Project (OpenRS).
21
22 OpenRS is free software: you can redistribute it and/or modify
23 it under the terms of the GNU General Public License as published by
24 the Free Software Foundation, either version 3 of the License, or
25 (at your option) any later version.
26
27 OpenRS is distributed in the hope that it will be useful,
28 but WITHOUT ANY WARRANTY; without even the implied warranty of
29 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30 GNU General Public License for more details.
31
32 You should have received a copy of the GNU General Public License
33 along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34*/
35
36#ifndef OPENRS_CFLCALCULATOR_HEADER
37#define OPENRS_CFLCALCULATOR_HEADER
38
39#include <opm/common/ErrorMacros.hpp>
40#include <opm/core/utility/Average.hpp>
41
44
45
46namespace Opm {
47 namespace cfl_calculator {
48
49
54 template <class Grid, class ReservoirProperties, class PressureSolution>
55 double findCFLtimeVelocity(const Grid& grid,
56 const ReservoirProperties& resprop,
57 const PressureSolution& pressure_sol)
58 {
59 double dt = 1e100;
60 typename Grid::CellIterator c = grid.cellbegin();
61 for (; c != grid.cellend(); ++c) {
62 double flux_p = 0.0;
63 double flux_n = 0.0;
64 typename Grid::CellIterator::FaceIterator f = c->facebegin();
65 for (; f != c->faceend(); ++f) {
66 const double loc_flux = pressure_sol.outflux(f);
67 if (loc_flux > 0) {
68 flux_p += loc_flux;
69 } else {
70 flux_n -= loc_flux;
71 }
72 }
73 double flux = std::max(flux_n, flux_p);
74 double loc_dt = (resprop.cflFactor()*c->volume()*resprop.porosity(c->index()))/flux;
75 if (loc_dt == 0.0) {
76 OPM_THROW(std::runtime_error, "Cfl computation gave dt = 0.0");
77 }
78 if (loc_dt < dt){
79 dt = loc_dt;
80 }
81 }
82 return dt;
83 }
84
85
90 template <class Grid, class ReservoirProperties>
91 double findCFLtimeGravity(const Grid& grid,
92 const ReservoirProperties& resprop,
93 const typename Grid::Vector& gravity)
94 {
95 typedef typename ReservoirProperties::PermTensor PermTensor;
96 typedef typename ReservoirProperties::MutablePermTensor MutablePermTensor;
97 const int dimension = Grid::Vector::dimension;
98 double dt = 1e100;
99 typename Grid::CellIterator c = grid.cellbegin();
100 for (; c != grid.cellend(); ++c) {
101 double flux = 0.0;
102 typename Grid::CellIterator::FaceIterator f = c->facebegin();
103 for (; f != c->faceend(); ++f) {
104 // UGLY WARNING
105 MutablePermTensor loc_perm_aver;
106 const double* permdata = 0;
107 if (!f->boundary()) {
108 PermTensor K0 = resprop.permeability(f->cellIndex());
109 PermTensor K1 = resprop.permeability(f->neighbourCellIndex());
110 loc_perm_aver = Opm::utils::arithmeticAverage<PermTensor, MutablePermTensor>(K0, K1);
111 permdata = loc_perm_aver.data();
112 } else {
113 permdata = resprop.permeability(f->cellIndex()).data();
114 }
115 PermTensor loc_perm(dimension, dimension, permdata);
116 typename Grid::Vector loc_halfface_normal = f->normal();
117 double loc_gravity_flux = 0.0;
118 for (int k = 0; k < dimension; ++k) {
119 for (int q = 0; q < dimension; ++q) {
120 loc_gravity_flux += loc_halfface_normal[q]*(loc_perm(q,k)*gravity[k]*resprop.densityDifference());
121 }
122 }
123 loc_gravity_flux *= f->area();
124 if (loc_gravity_flux > 0) {
125 flux += loc_gravity_flux;
126 }
127 }
128 double loc_dt = (resprop.cflFactorGravity()*c->volume()*resprop.porosity(c->index()))/flux;
129 if (loc_dt < dt){
130 dt = loc_dt;
131 }
132 }
133 return dt;
134 }
135
136
137
142 template <class Grid, class ReservoirProperties>
143 double findCFLtimeCapillary(const Grid& grid,
144 const ReservoirProperties& resprop)
145 {
146 typedef typename ReservoirProperties::PermTensor PermTensor;
147 typedef typename ReservoirProperties::MutablePermTensor MutablePermTensor;
148 const int dimension = Grid::Vector::dimension;
149 double dt = 1e100;
150 typename Grid::CellIterator c = grid.cellbegin();
151 for (; c != grid.cellend(); ++c) {
152 typename Grid::CellIterator::FaceIterator f = c->facebegin();
153 for (; f != c->faceend(); ++f) {
154 // UGLY WARNING
155 MutablePermTensor loc_perm_aver;
156 const double* permdata = 0;
157 if (!f->boundary()) {
158 PermTensor K0 = resprop.permeability(f->cellIndex());
159 PermTensor K1 = resprop.permeability(f->neighbourCellIndex());
160 loc_perm_aver = Opm::utils::arithmeticAverage<PermTensor, MutablePermTensor>(K0, K1);
161 permdata = loc_perm_aver.data();
162 } else {
163 permdata = resprop.permeability(f->cellIndex()).data();
164 }
165 // PermTensor loc_perm(dimension, dimension, permdata);
166 MutablePermTensor loc_perm(dimension, dimension, permdata);
167 MutablePermTensor loc_perm_inv = inverse3x3(loc_perm);
168 typename Grid::Vector loc_centroid = f->centroid();
169 loc_centroid -= c->centroid();
170 double spatial_contrib = loc_centroid*prod(loc_perm_inv, loc_centroid);
171 double loc_dt = spatial_contrib/resprop.cflFactorCapillary();
172 dt = std::min(dt, loc_dt);
173 }
174 }
175 return dt;
176 }
177
178 } // namespace cfl_calculator
179} // namespace Opm
180
181
182#endif // OPENRS_CFLCALCULATOR_HEADER
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
M inverse3x3(const M &m)
Definition: MatrixInverse.hpp:86
Dune::FieldVector< typename Matrix::value_type, rows > prod(const Matrix &A, const Dune::FieldVector< typename Matrix::value_type, rows > &x)
Matrix applied to a vector.
Definition: Matrix.hpp:669