mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Implemented VTK output for CpGrid using DUNE's VTKWriter and activated Matlab for CpGrid.
This commit is contained in:
parent
ec32758822
commit
c82778b3a9
@ -15,6 +15,9 @@
|
|||||||
|
|
||||||
#include <boost/filesystem.hpp>
|
#include <boost/filesystem.hpp>
|
||||||
|
|
||||||
|
#ifdef HAVE_DUNE_CORNERPOINT
|
||||||
|
#include <dune/grid/io/file/vtk/vtkwriter.hh>
|
||||||
|
#endif
|
||||||
namespace Opm
|
namespace Opm
|
||||||
{
|
{
|
||||||
|
|
||||||
@ -159,16 +162,35 @@ namespace Opm
|
|||||||
const int step,
|
const int step,
|
||||||
const std::string& output_dir)
|
const std::string& output_dir)
|
||||||
{
|
{
|
||||||
OPM_THROW(std::runtime_error, "outputStateVtk not implemented");
|
// Write data in VTK format.
|
||||||
}
|
std::ostringstream vtkfilename;
|
||||||
void outputStateMatlab(const Dune::CpGrid& grid,
|
std::ostringstream vtkpath;
|
||||||
const Opm::BlackoilState& state,
|
vtkpath << output_dir << "/vtk_files";
|
||||||
const int step,
|
vtkpath << "/output-" << std::setw(3) << std::setfill('0') << step;
|
||||||
const std::string& output_dir)
|
boost::filesystem::path fpath(vtkpath.str());
|
||||||
{
|
try {
|
||||||
OPM_THROW(std::runtime_error, "outputStateMatlab not implemented");
|
create_directories(fpath);
|
||||||
}
|
}
|
||||||
|
catch (...) {
|
||||||
|
OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath);
|
||||||
|
}
|
||||||
|
vtkfilename << "output-" << std::setw(3) << std::setfill('0') << step;
|
||||||
|
Dune::VTKWriter<Dune::CpGrid::LeafGridView> writer(grid.leafView(), Dune::VTK::nonconforming);
|
||||||
|
writer.addCellData(state.saturation(), "saturation", state.numPhases());
|
||||||
|
writer.addCellData(state.pressure(), "pressure", 1);
|
||||||
|
|
||||||
|
std::vector<double> cell_velocity;
|
||||||
|
Opm::estimateCellVelocity(AutoDiffGrid::numCells(grid),
|
||||||
|
AutoDiffGrid::numFaces(grid),
|
||||||
|
AutoDiffGrid::beginFaceCentroids(grid),
|
||||||
|
AutoDiffGrid::faceCells(grid),
|
||||||
|
AutoDiffGrid::beginCellCentroids(grid),
|
||||||
|
AutoDiffGrid::beginCellVolumes(grid),
|
||||||
|
AutoDiffGrid::dimensions(grid),
|
||||||
|
state.faceflux(), cell_velocity);
|
||||||
|
writer.addCellData(cell_velocity, "velocity", Dune::CpGrid::dimension);
|
||||||
|
writer.pwrite(vtkfilename.str(), vtkpath.str(), std::string("."), Dune::VTK::ascii);
|
||||||
|
}
|
||||||
#endif
|
#endif
|
||||||
|
|
||||||
}
|
}
|
||||||
|
@ -3,8 +3,17 @@
|
|||||||
#include <opm/core/grid.h>
|
#include <opm/core/grid.h>
|
||||||
#include <opm/core/simulator/BlackoilState.hpp>
|
#include <opm/core/simulator/BlackoilState.hpp>
|
||||||
#include <opm/core/simulator/WellState.hpp>
|
#include <opm/core/simulator/WellState.hpp>
|
||||||
|
#include <opm/core/utility/DataMap.hpp>
|
||||||
|
#include <opm/core/utility/miscUtilities.hpp>
|
||||||
|
|
||||||
|
#include <opm/autodiff/GridHelpers.hpp>
|
||||||
|
|
||||||
#include <string>
|
#include <string>
|
||||||
|
#include <sstream>
|
||||||
|
#include <iomanip>
|
||||||
|
#include <fstream>
|
||||||
|
|
||||||
|
#include <boost/filesystem.hpp>
|
||||||
|
|
||||||
#ifdef HAVE_DUNE_CORNERPOINT
|
#ifdef HAVE_DUNE_CORNERPOINT
|
||||||
#include <dune/grid/CpGrid.hpp>
|
#include <dune/grid/CpGrid.hpp>
|
||||||
@ -31,12 +40,51 @@ namespace Opm
|
|||||||
const Opm::BlackoilState& state,
|
const Opm::BlackoilState& state,
|
||||||
const int step,
|
const int step,
|
||||||
const std::string& output_dir);
|
const std::string& output_dir);
|
||||||
|
#endif
|
||||||
|
|
||||||
|
template<class Grid>
|
||||||
void outputStateMatlab(const Dune::CpGrid& grid,
|
void outputStateMatlab(const Grid& grid,
|
||||||
const Opm::BlackoilState& state,
|
const Opm::BlackoilState& state,
|
||||||
const int step,
|
const int step,
|
||||||
const std::string& output_dir);
|
const std::string& output_dir)
|
||||||
#endif
|
{
|
||||||
|
Opm::DataMap dm;
|
||||||
|
dm["saturation"] = &state.saturation();
|
||||||
|
dm["pressure"] = &state.pressure();
|
||||||
|
dm["surfvolume"] = &state.surfacevol();
|
||||||
|
std::vector<double> cell_velocity;
|
||||||
|
Opm::estimateCellVelocity(AutoDiffGrid::numCells(grid),
|
||||||
|
AutoDiffGrid::numFaces(grid),
|
||||||
|
AutoDiffGrid::beginFaceCentroids(grid),
|
||||||
|
AutoDiffGrid::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
|
#endif
|
||||||
|
Loading…
Reference in New Issue
Block a user