36#ifndef OPENRS_SIMULATORUTILITIES_HEADER
37#define OPENRS_SIMULATORUTILITIES_HEADER
40#include <opm/common/utility/platform_dependent/disable_warnings.h>
42#include <dune/common/version.hh>
43#include <dune/common/fvector.hh>
44#include <dune/grid/io/file/vtk/vtkwriter.hh>
46#include <opm/common/utility/platform_dependent/reenable_warnings.h>
48#include <opm/common/ErrorMacros.hpp>
64 template <
class Gr
idInterface,
class FlowSol>
67 const FlowSol& flow_solution)
71 cell_velocity.clear();
72 cell_velocity.resize(ginterf.numberOfCells());
73 for (
typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
75 typename GridInterface::Vector cell_v(0.0);
76 typename GridInterface::CellIterator::FaceIterator f = c->facebegin();
77 for (; f != c->faceend(); ++f, ++numf) {
78 double flux = flow_solution.outflux(f);
79 typename GridInterface::Vector v = f->centroid();
81 v *= flux/c->volume();
84 cell_velocity[c->index()] = cell_v;
92 template <
class Gr
idInterface>
95 const std::vector<double>& face_flux)
99 typedef typename GridInterface::Vector Vec;
100 cell_velocity.clear();
101 cell_velocity.resize(grid.numCells(), Vec(0.0));
102 for (
int face = 0; face < grid.numFaces(); ++face) {
103 int c[2] = { grid.faceCell(face, 0), grid.faceCell(face, 1) };
104 Vec fc = grid.faceCentroid(face);
105 double flux = face_flux[face];
106 for (
int i = 0; i < 2; ++i) {
108 Vec v_contrib = fc - grid.cellCentroid(c[i]);
109 v_contrib *= flux/grid.cellVolume(c[i]);
110 cell_velocity[c[i]] += (i == 0) ? v_contrib : -v_contrib;
125 template <
class Gr
idInterface,
class FlowSol>
128 const FlowSol& flow_solution,
129 const std::vector<int>& partition,
130 const int my_partition)
134 cell_velocity.clear();
135 cell_velocity.resize(ginterf.numberOfCells());
136 for (
typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
137 if (partition[c->index()] != my_partition) {
138 cell_velocity[c->index()] = 0.0;
141 typename GridInterface::Vector cell_v(0.0);
142 typename GridInterface::CellIterator::FaceIterator f = c->facebegin();
143 for (; f != c->faceend(); ++f, ++numf) {
144 double flux = flow_solution.outflux(f);
145 typename GridInterface::Vector v = f->centroid();
147 v *= flux/c->volume();
150 cell_velocity[c->index()] = cell_v;
156 template <
class ReservoirProperty>
157 void computePhaseVelocities(std::vector<Dune::FieldVector<double, 3> >& phase_velocity_water,
158 std::vector<Dune::FieldVector<double, 3> >& phase_velocity_oil,
159 const ReservoirProperty& res_prop,
160 const std::vector<double>& saturation,
161 const std::vector<Dune::FieldVector<double, 3> >& cell_velocity)
163 assert(saturation.size() == cell_velocity.size());
164 int num_cells = saturation.size();
165 phase_velocity_water = cell_velocity;
166 phase_velocity_oil = cell_velocity;
167 for (
int i = 0; i < num_cells; ++i) {
168 double f = res_prop.fractionalFlow(i, saturation[i]);
169 phase_velocity_water[i] *= f;
170 phase_velocity_oil[i] *= (1.0 - f);
181 template <
class Gr
idInterface,
class FlowSol>
184 const FlowSol& flow_solution)
186 cell_pressure.clear();
187 cell_pressure.resize(ginterf.numberOfCells());
188 for (
typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
189 cell_pressure[c->index()] = flow_solution.pressure(c);
197 template <
class Gr
idInterface,
class FlowSol>
200 const FlowSol& flow_solution,
201 const std::vector<int>& partition,
202 const int my_partition)
204 cell_pressure.clear();
205 cell_pressure.resize(ginterf.numberOfCells());
206 for (
typename GridInterface::CellIterator c = ginterf.cellbegin(); c != ginterf.cellend(); ++c) {
207 if (partition[c->index()] != my_partition) {
208 cell_pressure[c->index()] = 0.0;
210 cell_pressure[c->index()] = flow_solution.pressure(c);
222 template <
class ReservoirProperties>
224 const ReservoirProperties& rp,
225 const std::vector<double>& sat)
227 int num_cells = sat.size();
228 cap_pressure.resize(num_cells);
229 for (
int cell = 0; cell < num_cells; ++cell) {
230 cap_pressure[cell] = rp.capillaryPressure(cell, sat[cell]);
236 template <
class Gr
idInterface,
class ReservoirProperties,
class FlowSol>
237 void writeVtkOutput(
const GridInterface& ginterf,
238 const ReservoirProperties& rp,
239 const FlowSol& flowsol,
240 const std::vector<double>& saturation,
241 const std::string& filename)
244 typedef typename GridInterface::Vector Vec;
245 std::vector<Vec> cell_velocity;
247 std::array<std::vector<Vec>, 2> phase_velocities;
248 computePhaseVelocities(phase_velocities[0], phase_velocities[1], rp, saturation, cell_velocity);
250 std::vector<double> cell_velocity_flat(&*cell_velocity.front().begin(),
251 &*cell_velocity.back().end());
252 std::vector<double> water_velocity_flat(&*phase_velocities[0].front().begin(),
253 &*phase_velocities[0].back().end());
254 std::vector<double> oil_velocity_flat(&*phase_velocities[1].front().begin(),
255 &*phase_velocities[1].back().end());
256 std::vector<double> cell_pressure;
258 std::vector<double> cap_pressure;
260 int num_cells = saturation.size();
264 std::vector<double> fractional_flow_(num_cells);
265 for (
int i = 0; i < num_cells; ++i) {
269 fractional_flow_[i] = rp.fractionalFlow(i, saturation[i]);
273 Dune::VTKWriter<typename GridInterface::GridType::LeafGridView> vtkwriter(ginterf.grid().leafGridView());
274 vtkwriter.addCellData(saturation,
"saturation");
275 vtkwriter.addCellData(cell_pressure,
"pressure");
276 vtkwriter.addCellData(cap_pressure,
"capillary pressure");
277 vtkwriter.addCellData(fractional_flow_,
"fractional flow [water]");
280 vtkwriter.addCellData(cell_velocity_flat,
"velocity", Vec::dimension);
281 vtkwriter.addCellData(water_velocity_flat,
"phase velocity [water]", Vec::dimension);
282 vtkwriter.addCellData(oil_velocity_flat,
"phase velocity [oil]", Vec::dimension);
283 vtkwriter.write(filename, Dune::VTK::ascii);
287 inline void writeField(
const std::vector<double>& field,
288 const std::string& filename)
290 std::ofstream os(filename.c_str());
292 OPM_THROW(std::runtime_error,
"Could not open file " + filename);
294 os << field.size() <<
'\n';
295 std::ranges::copy(field, std::ostream_iterator<double>(os,
"\n"));
Inverting small matrices.
Definition ImplicitAssembly.hpp:43
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:93
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:223
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:65
void getCellPressure(std::vector< double > &cell_pressure, const GridInterface &ginterf, const FlowSol &flow_solution)
Definition SimulatorUtilities.hpp:182