SimulatorUtilities.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: SimulatorUtilities.hpp
4//
5// Created: Fri Aug 28 15:00:15 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Bård Skaflestad <bard.skaflestad@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_SIMULATORUTILITIES_HEADER
37#define OPENRS_SIMULATORUTILITIES_HEADER
38
39
40#include <opm/common/ErrorMacros.hpp>
41#include <dune/common/version.hh>
42#include <dune/common/fvector.hh>
43#include <dune/grid/io/file/vtk/vtkwriter.hh>
44#include <vector>
45#include <fstream>
46#include <algorithm>
47#include <iterator>
48
49namespace Opm
50{
51
52
59 template <class GridInterface, class FlowSol>
60 void estimateCellVelocity(std::vector<typename GridInterface::Vector>& cell_velocity,
61 const GridInterface& ginterf,
62 const FlowSol& flow_solution)
63 {
64 // Algorithm used is same as in halfFaceFluxToCellVelocity.hpp
65 // in the Sintef legacy c++ code.
66 cell_velocity.clear();
67 cell_velocity.resize(ginterf.numberOfCells());
68 for (typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
69 int numf = 0;
70 typename GridInterface::Vector cell_v(0.0);
71 typename GridInterface::CellIterator::FaceIterator f = c->facebegin();
72 for (; f != c->faceend(); ++f, ++numf) {
73 double flux = flow_solution.outflux(f);
74 typename GridInterface::Vector v = f->centroid();
75 v -= c->centroid();
76 v *= flux/c->volume();
77 cell_v += v;
78 }
79 cell_velocity[c->index()] = cell_v;//.two_norm();
80 }
81 }
82
89 template <class GridInterface>
90 void estimateCellVelocitySimpleInterface(std::vector<typename GridInterface::Vector>& cell_velocity,
91 const GridInterface& grid,
92 const std::vector<double>& face_flux)
93 {
94 // Algorithm used is same as in halfFaceFluxToCellVelocity.hpp
95 // in the Sintef legacy c++ code.
96 typedef typename GridInterface::Vector Vec;
97 cell_velocity.clear();
98 cell_velocity.resize(grid.numCells(), Vec(0.0));
99 for (int face = 0; face < grid.numFaces(); ++face) {
100 int c[2] = { grid.faceCell(face, 0), grid.faceCell(face, 1) };
101 Vec fc = grid.faceCentroid(face);
102 double flux = face_flux[face];
103 for (int i = 0; i < 2; ++i) {
104 if (c[i] >= 0) {
105 Vec v_contrib = fc - grid.cellCentroid(c[i]);
106 v_contrib *= flux/grid.cellVolume(c[i]);
107 cell_velocity[c[i]] += (i == 0) ? v_contrib : -v_contrib;
108 }
109 }
110 }
111 }
112
113
122 template <class GridInterface, class FlowSol>
123 void estimateCellVelocity(std::vector<typename GridInterface::Vector>& cell_velocity,
124 const GridInterface& ginterf,
125 const FlowSol& flow_solution,
126 const std::vector<int>& partition,
127 const int my_partition)
128 {
129 // Algorithm used is same as in halfFaceFluxToCellVelocity.hpp
130 // in the Sintef legacy c++ code.
131 cell_velocity.clear();
132 cell_velocity.resize(ginterf.numberOfCells());
133 for (typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
134 if (partition[c->index()] != my_partition) {
135 cell_velocity[c->index()] = 0.0;
136 } else {
137 int numf = 0;
138 typename GridInterface::Vector cell_v(0.0);
139 typename GridInterface::CellIterator::FaceIterator f = c->facebegin();
140 for (; f != c->faceend(); ++f, ++numf) {
141 double flux = flow_solution.outflux(f);
142 typename GridInterface::Vector v = f->centroid();
143 v -= c->centroid();
144 v *= flux/c->volume();
145 cell_v += v;
146 }
147 cell_velocity[c->index()] = cell_v;//.two_norm();
148 }
149 }
150 }
151
152
153 template <class ReservoirProperty>
154 void computePhaseVelocities(std::vector<Dune::FieldVector<double, 3> >& phase_velocity_water,
155 std::vector<Dune::FieldVector<double, 3> >& phase_velocity_oil,
156 const ReservoirProperty& res_prop,
157 const std::vector<double>& saturation,
158 const std::vector<Dune::FieldVector<double, 3> >& cell_velocity)
159 {
160 assert(saturation.size() == cell_velocity.size());
161 int num_cells = saturation.size();
162 phase_velocity_water = cell_velocity;
163 phase_velocity_oil = cell_velocity;
164 for (int i = 0; i < num_cells; ++i) {
165 double f = res_prop.fractionalFlow(i, saturation[i]);
166 phase_velocity_water[i] *= f;
167 phase_velocity_oil[i] *= (1.0 - f);
168 }
169 }
170
171
172
173
178 template <class GridInterface, class FlowSol>
179 void getCellPressure(std::vector<double>& cell_pressure,
180 const GridInterface& ginterf,
181 const FlowSol& flow_solution)
182 {
183 cell_pressure.clear();
184 cell_pressure.resize(ginterf.numberOfCells());
185 for (typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
186 cell_pressure[c->index()] = flow_solution.pressure(c);
187 }
188 }
189
194 template <class GridInterface, class FlowSol>
195 void getCellPressure(std::vector<double>& cell_pressure,
196 const GridInterface& ginterf,
197 const FlowSol& flow_solution,
198 const std::vector<int>& partition,
199 const int my_partition)
200 {
201 cell_pressure.clear();
202 cell_pressure.resize(ginterf.numberOfCells());
203 for (typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
204 if (partition[c->index()] != my_partition) {
205 cell_pressure[c->index()] = 0.0;
206 } else {
207 cell_pressure[c->index()] = flow_solution.pressure(c);
208 }
209 }
210 }
211
212
213
219 template <class ReservoirProperties>
220 void computeCapPressure(std::vector<double>& cap_pressure,
221 const ReservoirProperties& rp,
222 const std::vector<double>& sat)
223 {
224 int num_cells = sat.size();
225 cap_pressure.resize(num_cells);
226 for (int cell = 0; cell < num_cells; ++cell) {
227 cap_pressure[cell] = rp.capillaryPressure(cell, sat[cell]);
228 }
229 }
230
231
233 template <class GridInterface, class ReservoirProperties, class FlowSol>
234 void writeVtkOutput(const GridInterface& ginterf,
235 const ReservoirProperties& rp,
236 const FlowSol& flowsol,
237 const std::vector<double>& saturation,
238 const std::string& filename)
239 {
240 // Extract data in proper format.
241 typedef typename GridInterface::Vector Vec;
242 std::vector<Vec> cell_velocity;
243 estimateCellVelocity(cell_velocity, ginterf, flowsol);
244 std::array<std::vector<Vec>, 2> phase_velocities;
245 computePhaseVelocities(phase_velocities[0], phase_velocities[1], rp, saturation, cell_velocity);
246 // Dune's vtk writer wants multi-component data to be flattened.
247 std::vector<double> cell_velocity_flat(&*cell_velocity.front().begin(),
248 &*cell_velocity.back().end());
249 std::vector<double> water_velocity_flat(&*phase_velocities[0].front().begin(),
250 &*phase_velocities[0].back().end());
251 std::vector<double> oil_velocity_flat(&*phase_velocities[1].front().begin(),
252 &*phase_velocities[1].back().end());
253 std::vector<double> cell_pressure;
254 getCellPressure(cell_pressure, ginterf, flowsol);
255 std::vector<double> cap_pressure;
256 computeCapPressure(cap_pressure, rp, saturation);
257 int num_cells = saturation.size();
258// std::array<std::vector<double>, 2> phase_mobilities_;
259// phase_mobilities_[0].resize(num_cells);
260// phase_mobilities_[1].resize(num_cells);
261 std::vector<double> fractional_flow_(num_cells);
262 for (int i = 0; i < num_cells; ++i) {
263// for (int phase = 0; phase < 2; ++phase) {
264// rp.phaseMobility(phase, i, saturation[i], phase_mobilities_[phase][i]);
265// }
266 fractional_flow_[i] = rp.fractionalFlow(i, saturation[i]);
267 }
268
269 // Write data.
270#if DUNE_VERSION_NEWER(DUNE_GRID, 2, 3)
271 Dune::VTKWriter<typename GridInterface::GridType::LeafGridView> vtkwriter(ginterf.grid().leafGridView());
272#else
273 Dune::VTKWriter<typename GridInterface::GridType::LeafGridView> vtkwriter(ginterf.grid().leafView());
274#endif
275 vtkwriter.addCellData(saturation, "saturation");
276 vtkwriter.addCellData(cell_pressure, "pressure");
277 vtkwriter.addCellData(cap_pressure, "capillary pressure");
278 vtkwriter.addCellData(fractional_flow_, "fractional flow [water]");
279// vtkwriter.addCellData(phase_mobilities_[0], "phase mobility [water]");
280// vtkwriter.addCellData(phase_mobilities_[1], "phase mobility [oil]");
281 vtkwriter.addCellData(cell_velocity_flat, "velocity", Vec::dimension);
282 vtkwriter.addCellData(water_velocity_flat, "phase velocity [water]", Vec::dimension);
283 vtkwriter.addCellData(oil_velocity_flat, "phase velocity [oil]", Vec::dimension);
284 vtkwriter.write(filename, Dune::VTK::ascii);
285 }
286
287
288 inline void writeField(const std::vector<double>& field,
289 const std::string& filename)
290 {
291 std::ofstream os(filename.c_str());
292 if (!os) {
293 OPM_THROW(std::runtime_error, "Could not open file " << filename);
294 }
295 os << field.size() << '\n';
296 std::copy(field.begin(), field.end(), std::ostream_iterator<double>(os, "\n"));
297 }
298
299} // namespace Opm
300
301
302#endif // OPENRS_SIMULATORUTILITIES_HEADER
Definition: BlackoilFluid.hpp:32
void computePhaseVelocities(std::vector< Dune::FieldVector< double, 3 > > &phase_velocity_water, std::vector< Dune::FieldVector< double, 3 > > &phase_velocity_oil, const ReservoirProperty &res_prop, const std::vector< double > &saturation, const std::vector< Dune::FieldVector< double, 3 > > &cell_velocity)
Definition: SimulatorUtilities.hpp:154
void writeField(const std::vector< double > &field, const std::string &filename)
Definition: SimulatorUtilities.hpp:288
void writeVtkOutput(const GridInterface &ginterf, const ReservoirProperties &rp, const FlowSol &flowsol, const std::vector< double > &saturation, const std::string &filename)
Definition: SimulatorUtilities.hpp:234
void estimateCellVelocitySimpleInterface(std::vector< typename GridInterface::Vector > &cell_velocity, const GridInterface &grid, const std::vector< double > &face_flux)
Estimates a scalar cell velocity from face fluxes.
Definition: SimulatorUtilities.hpp:90
void computeCapPressure(std::vector< double > &cap_pressure, const ReservoirProperties &rp, const std::vector< double > &sat)
Computes the capillary pressure in each cell from the cell saturations.
Definition: SimulatorUtilities.hpp:220
void estimateCellVelocity(std::vector< typename GridInterface::Vector > &cell_velocity, const GridInterface &ginterf, const FlowSol &flow_solution)
Estimates a scalar cell velocity from outgoing fluxes.
Definition: SimulatorUtilities.hpp:60
void getCellPressure(std::vector< double > &cell_pressure, const GridInterface &ginterf, const FlowSol &flow_solution)
Definition: SimulatorUtilities.hpp:179