Add "output_vtk" parameter, split output function in two.

Vtk and Matlab output now happens in two different functions.
This commit is contained in:
Atgeirr Flø Rasmussen 2012-09-17 07:54:50 +02:00
parent 3d883a5aed
commit 076164d3f3

View File

@ -64,7 +64,11 @@ namespace Opm
namespace
{
void outputState(const UnstructuredGrid& grid,
void outputStateVtk(const UnstructuredGrid& grid,
const Opm::PolymerState& state,
const int step,
const std::string& output_dir);
void outputStateMatlab(const UnstructuredGrid& grid,
const Opm::PolymerState& state,
const int step,
const std::string& output_dir);
@ -72,7 +76,6 @@ namespace Opm
const std::string& output_dir);
void outputWellReport(const Opm::WellReport& wellreport,
const std::string& output_dir);
} // anonymous namespace
@ -100,6 +103,7 @@ namespace Opm
// Parameters for output.
bool output_;
bool output_vtk_;
std::string output_dir_;
int output_interval_;
// Parameters for transport solver.
@ -191,6 +195,7 @@ namespace Opm
// For output.
output_ = param.getDefault("output", true);
if (output_) {
output_vtk_ = param.getDefault("output_vtk", true);
output_dir_ = param.getDefault("output_dir", std::string("output"));
// Ensure that output dir exists
boost::filesystem::path fpath(output_dir_);
@ -283,7 +288,10 @@ namespace Opm
// Report timestep and (optionally) write state to disk.
timer.report(std::cout);
if (output_ && (timer.currentStepNum() % output_interval_ == 0)) {
outputState(grid_, state, timer.currentStepNum(), output_dir_);
if (output_vtk_) {
outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
}
outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
}
// Solve pressure.
@ -398,7 +406,10 @@ namespace Opm
}
if (output_) {
outputState(grid_, state, timer.currentStepNum(), output_dir_);
if (output_vtk_) {
outputStateVtk(grid_, state, timer.currentStepNum(), output_dir_);
}
outputStateMatlab(grid_, state, timer.currentStepNum(), output_dir_);
outputWaterCut(watercut, output_dir_);
if (wells_) {
outputWellReport(wellreport, output_dir_);
@ -423,14 +434,22 @@ namespace Opm
namespace
{
void outputState(const UnstructuredGrid& grid,
void outputStateVtk(const UnstructuredGrid& grid,
const Opm::PolymerState& state,
const int step,
const std::string& output_dir)
{
// Write data in VTK format.
std::ostringstream vtkfilename;
vtkfilename << output_dir << "/output-" << std::setw(3) << std::setfill('0') << step << ".vtu";
vtkfilename << output_dir << "/vtk_files";
boost::filesystem::path fpath(vtkfilename.str());
try {
create_directories(fpath);
}
catch (...) {
THROW("Creating directories failed: " << fpath);
}
vtkfilename << "/output-" << std::setw(3) << std::setfill('0') << step << ".vtu";
std::ofstream vtkfile(vtkfilename.str().c_str());
if (!vtkfile) {
THROW("Failed to open " << vtkfilename.str());
@ -458,6 +477,40 @@ namespace Opm
}
}
void outputStateMatlab(const UnstructuredGrid& grid,
const Opm::PolymerState& state,
const int step,
const std::string& output_dir)
{
Opm::DataMap dm;
dm["saturation"] = &state.saturation();
dm["pressure"] = &state.pressure();
dm["concentration"] = &state.concentration();
dm["cmax"] = &state.maxconcentration();
std::vector<double> cell_velocity;
Opm::estimateCellVelocity(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 (...) {
THROW("Creating directories failed: " << fpath);
}
fname << "/" << std::setw(3) << std::setfill('0') << step << ".txt";
std::ofstream file(fname.str().c_str());
if (!file) {
THROW("Failed to open " << fname.str());
}
const std::vector<double>& d = *(it->second);
std::copy(d.begin(), d.end(), std::ostream_iterator<double>(file, "\n"));
}
}
void outputWaterCut(const Opm::Watercut& watercut,
const std::string& output_dir)