opm-simulators/opm/autodiff/SimulatorFullyImplicitBlackoilOutput.hpp
Markus Blatt 7127101c1c Makes distinction between functions more clear.
Currently, there are two abstract interface for the grids. One that
usually returns pods and arrays of them that also can be used by C
and is used also in opm-core, and one that returns Eigen datastructures
 needed within opm-autodiff.

This commit adds a postfix ToEigen to those functions (faceCells, and
cellCentroidsZ) one could imagine to also return pods and arrays of them.
This should at least resolve the confusion about the two faceCells functions.

The next step will be issue #192
Fixes #176
2014-08-28 14:44:13 +02:00

95 lines
3.4 KiB
C++

#ifndef OPM_SIMULATORFULLYIMPLICITBLACKOILOUTPUT_HEADER_INCLUDED
#define OPM_SIMULATORFULLYIMPLICITBLACKOILOUTPUT_HEADER_INCLUDED
#include <opm/core/grid.h>
#include <opm/core/simulator/BlackoilState.hpp>
#include <opm/core/simulator/WellState.hpp>
#include <opm/core/utility/DataMap.hpp>
#include <opm/core/utility/ErrorMacros.hpp>
#include <opm/core/utility/miscUtilities.hpp>
#include <opm/autodiff/GridHelpers.hpp>
#include <string>
#include <sstream>
#include <iomanip>
#include <fstream>
#include <boost/filesystem.hpp>
#ifdef HAVE_DUNE_CORNERPOINT
#include <dune/grid/CpGrid.hpp>
#endif
namespace Opm
{
void outputStateVtk(const UnstructuredGrid& grid,
const Opm::BlackoilState& state,
const int step,
const std::string& output_dir);
void outputStateMatlab(const UnstructuredGrid& grid,
const Opm::BlackoilState& state,
const int step,
const std::string& output_dir);
void outputWellStateMatlab(const Opm::WellState& well_state,
const int step,
const std::string& output_dir);
#ifdef HAVE_DUNE_CORNERPOINT
void outputStateVtk(const Dune::CpGrid& grid,
const Opm::BlackoilState& state,
const int step,
const std::string& output_dir);
#endif
template<class Grid>
void outputStateMatlab(const Grid& grid,
const Opm::BlackoilState& state,
const int step,
const std::string& output_dir)
{
Opm::DataMap dm;
dm["saturation"] = &state.saturation();
dm["pressure"] = &state.pressure();
dm["surfvolume"] = &state.surfacevol();
dm["rs"] = &state.gasoilratio();
dm["rv"] = &state.rv();
std::vector<double> cell_velocity;
Opm::estimateCellVelocity(AutoDiffGrid::numCells(grid),
AutoDiffGrid::numFaces(grid),
AutoDiffGrid::beginFaceCentroids(grid),
UgGridHelpers::faceCells(grid),
AutoDiffGrid::beginCellCentroids(grid),
AutoDiffGrid::beginCellVolumes(grid),
AutoDiffGrid::dimensions(grid),
state.faceflux(), cell_velocity);
dm["velocity"] = &cell_velocity;
// Write data (not grid) in Matlab format
for (Opm::DataMap::const_iterator it = dm.begin(); it != dm.end(); ++it) {
std::ostringstream fname;
fname << output_dir << "/" << it->first;
boost::filesystem::path fpath = fname.str();
try {
create_directories(fpath);
}
catch (...) {
OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath);
}
fname << "/" << std::setw(3) << std::setfill('0') << step << ".txt";
std::ofstream file(fname.str().c_str());
if (!file) {
OPM_THROW(std::runtime_error, "Failed to open " << fname.str());
}
file.precision(15);
const std::vector<double>& d = *(it->second);
std::copy(d.begin(), d.end(), std::ostream_iterator<double>(file, "\n"));
}
}
}
#endif