/* Copyright (c) 2014 SINTEF ICT, Applied Mathematics. Copyright (c) 2015-2016 IRIS AS This file is part of the Open Porous Media project (OPM). OPM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OPM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OPM. If not, see . */ #include "config.h" #include "SimulatorFullyImplicitBlackoilOutput.hpp" #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include //For OutputWriterHelper #include #include #ifdef HAVE_OPM_GRID #include #include #include #include #endif namespace Opm { void outputStateVtk(const UnstructuredGrid& grid, const SimulationDataContainer& state, const int step, const std::string& output_dir) { // Write data in VTK format. std::ostringstream vtkfilename; vtkfilename << output_dir << "/vtk_files"; boost::filesystem::path fpath(vtkfilename.str()); try { create_directories(fpath); } catch (...) { OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath); } vtkfilename << "/output-" << std::setw(3) << std::setfill('0') << step << ".vtu"; std::ofstream vtkfile(vtkfilename.str().c_str()); if (!vtkfile) { OPM_THROW(std::runtime_error, "Failed to open " << vtkfilename.str()); } Opm::DataMap dm; dm["saturation"] = &state.saturation(); dm["pressure"] = &state.pressure(); std::vector 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; Opm::writeVtkData(grid, dm, vtkfile); } void outputStateMatlab(const UnstructuredGrid& grid, const Opm::SimulationDataContainer& state, const int step, const std::string& output_dir) { Opm::DataMap dm; dm["saturation"] = &state.saturation(); dm["pressure"] = &state.pressure(); for (const auto& pair : state.cellData()) { const std::string& name = pair.first; std::string key; if( name == "SURFACEVOL" ) { key = "surfvolume"; } else if( name == "RV" ) { key = "rv"; } else if( name == "GASOILRATIO" ) { key = "rs"; } else { // otherwise skip entry continue; } // set data to datmap dm[ key ] = &pair.second; } std::vector 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& d = *(it->second); std::copy(d.begin(), d.end(), std::ostream_iterator(file, "\n")); } } void outputWellStateMatlab(const Opm::WellState& well_state, const int step, const std::string& output_dir) { Opm::DataMap dm; dm["bhp"] = &well_state.bhp(); dm["wellrates"] = &well_state.wellRates(); // 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& d = *(it->second); std::copy(d.begin(), d.end(), std::ostream_iterator(file, "\n")); } } #if 0 void outputWaterCut(const Opm::Watercut& watercut, const std::string& output_dir) { // Write water cut curve. std::string fname = output_dir + "/watercut.txt"; std::ofstream os(fname.c_str()); if (!os) { OPM_THROW(std::runtime_error, "Failed to open " << fname); } watercut.write(os); } void outputWellReport(const Opm::WellReport& wellreport, const std::string& output_dir) { // Write well report. std::string fname = output_dir + "/wellreport.txt"; std::ofstream os(fname.c_str()); if (!os) { OPM_THROW(std::runtime_error, "Failed to open " << fname); } wellreport.write(os); } #endif #ifdef HAVE_OPM_GRID void outputStateVtk(const Dune::CpGrid& grid, const Opm::SimulationDataContainer& state, const int step, const std::string& output_dir) { // Write data in VTK format. std::ostringstream vtkfilename; std::ostringstream vtkpath; vtkpath << output_dir << "/vtk_files"; vtkpath << "/output-" << std::setw(3) << std::setfill('0') << step; boost::filesystem::path fpath(vtkpath.str()); try { create_directories(fpath); } catch (...) { OPM_THROW(std::runtime_error, "Creating directories failed: " << fpath); } vtkfilename << "output-" << std::setw(3) << std::setfill('0') << step; #if DUNE_VERSION_NEWER(DUNE_GRID, 2, 3) Dune::VTKWriter writer(grid.leafGridView(), Dune::VTK::nonconforming); #else Dune::VTKWriter writer(grid.leafView(), Dune::VTK::nonconforming); #endif writer.addCellData(state.saturation(), "saturation", state.numPhases()); writer.addCellData(state.pressure(), "pressure", 1); std::vector 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 namespace detail { struct WriterCall : public ThreadHandle :: ObjectInterface { BlackoilOutputWriter& writer_; std::unique_ptr< SimulatorTimerInterface > timer_; const SimulationDataContainer state_; const WellState wellState_; std::vector simProps_; const bool substep_; explicit WriterCall( BlackoilOutputWriter& writer, const SimulatorTimerInterface& timer, const SimulationDataContainer& state, const WellState& wellState, const std::vector& simProps, bool substep ) : writer_( writer ), timer_( timer.clone() ), state_( state ), wellState_( wellState ), simProps_( simProps ), substep_( substep ) { } // callback to writer's serial writeTimeStep method void run () { // write data writer_.writeTimeStepSerial( *timer_, state_, wellState_, simProps_, substep_ ); } }; } void BlackoilOutputWriter:: writeTimeStepWithoutCellProperties( const SimulatorTimerInterface& timer, const SimulationDataContainer& localState, const WellState& localWellState, bool substep) { std::vector noCellProperties; writeTimeStepWithCellProperties(timer, localState, localWellState, noCellProperties, substep); } void BlackoilOutputWriter:: writeTimeStepWithCellProperties( const SimulatorTimerInterface& timer, const SimulationDataContainer& localState, const WellState& localWellState, const std::vector& cellData, bool substep) { // VTK output (is parallel if grid is parallel) if( vtkWriter_ ) { vtkWriter_->writeTimeStep( timer, localState, localWellState, false ); } bool isIORank = output_ ; if( parallelOutput_ && parallelOutput_->isParallel() ) { // If this is not the initial write and no substep, then the well // state used in the computation is actually the one of the last // step. We need that well state for the gathering. Otherwise // It an exception with a message like "global state does not // contain well ..." might be thrown. int wellStateStepNumber = ( ! substep && timer.reportStepNum() > 0) ? (timer.reportStepNum() - 1) : timer.reportStepNum(); // collect all solutions to I/O rank isIORank = parallelOutput_->collectToIORank( localState, localWellState, wellStateStepNumber ); } const SimulationDataContainer& state = (parallelOutput_ && parallelOutput_->isParallel() ) ? parallelOutput_->globalReservoirState() : localState; const WellState& wellState = (parallelOutput_ && parallelOutput_->isParallel() ) ? parallelOutput_->globalWellState() : localWellState; // serial output is only done on I/O rank if( isIORank ) { if( asyncOutput_ ) { // dispatch the write call to the extra thread asyncOutput_->dispatch( detail::WriterCall( *this, timer, state, wellState, cellData, substep ) ); } else { // just write the data to disk writeTimeStepSerial( timer, state, wellState, cellData, substep ); } } } void BlackoilOutputWriter:: writeTimeStepSerial(const SimulatorTimerInterface& timer, const SimulationDataContainer& state, const WellState& wellState, const std::vector& simProps, bool substep) { // Matlab output if( matlabWriter_ ) { matlabWriter_->writeTimeStep( timer, state, wellState, substep ); } // ECL output if ( eclWriter_ ) { const auto& initConfig = eclipseState_->getInitConfig(); if (initConfig.restartRequested() && ((initConfig.getRestartStep()) == (timer.currentStepNum()))) { std::cout << "Skipping restart write in start of step " << timer.currentStepNum() << std::endl; } else { eclWriter_->writeTimeStep(timer.reportStepNum(), substep, timer.simulationTimeElapsed(), simToSolution( state, phaseUsage_ ), wellState.report(phaseUsage_), simProps); } } // write backup file if( backupfile_.is_open() ) { int reportStep = timer.reportStepNum(); int currentTimeStep = timer.currentStepNum(); if( (reportStep == currentTimeStep || // true for SimulatorTimer currentTimeStep == 0 || // true for AdaptiveSimulatorTimer at reportStep timer.done() ) // true for AdaptiveSimulatorTimer at reportStep && lastBackupReportStep_ != reportStep ) // only backup report step once { // store report step lastBackupReportStep_ = reportStep; // write resport step number backupfile_.write( (const char *) &reportStep, sizeof(int) ); try { backupfile_ << state; const WellStateFullyImplicitBlackoil& boWellState = static_cast< const WellStateFullyImplicitBlackoil& > (wellState); backupfile_ << boWellState; } catch ( const std::bad_cast& e ) { } backupfile_ << std::flush; } } // end backup } void BlackoilOutputWriter:: restore(SimulatorTimerInterface& timer, BlackoilState& state, WellStateFullyImplicitBlackoil& wellState, const std::string& filename, const int desiredResportStep ) { std::ifstream restorefile( filename.c_str() ); if( restorefile ) { std::cout << "============================================================================"<getInitConfig(); return initconfig.restartRequested(); } }