// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- // vi: set et ts=4 sw=4 sts=4: /* 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 2 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 . Consult the COPYING file in the top-level source directory of this module for the precise wording of the license and the list of copyright holders. */ /*! * \file * * \copydoc Opm::EclWriter */ #ifndef EWOMS_ECL_GENERIC_WRITER_HH #define EWOMS_ECL_GENERIC_WRITER_HH #include "collecttoiorank.hh" #include #if HAVE_DUNE_ALUGRID #include "alucartesianindexmapper.hh" #include #endif // HAVE_DUNE_ALUGRID #include #include #include #include #include #include #include namespace Opm { class EclipseIO; class EclipseState; class EclInterRegFlowMap; class Inplace; struct NNCdata; class Schedule; class SummaryConfig; class SummaryState; class UDQState; } // namespace Opm namespace Opm { namespace Action { class State; }} // namespace Opm::Action namespace Opm { template class EclGenericWriter { typedef Dune::CartesianIndexMapper CartesianIndexMapper; typedef Dune::CartesianIndexMapper EquilCartesianIndexMapper; using CollectDataToIORankType = CollectDataToIORank; using TransmissibilityType = EclTransmissibility; public: // The Simulator object should preferably have been const - the // only reason that is not the case is due to the SummaryState // object owned deep down by the vanguard. EclGenericWriter(const Schedule& schedule, const EclipseState& eclState, const SummaryConfig& summaryConfig, const Grid& grid, const EquilGrid* equilGrid, const GridView& gridView, const Dune::CartesianIndexMapper& cartMapper, const Dune::CartesianIndexMapper* equilCartMapper, bool enableAsyncOutput, bool enableEsmry); const EclipseIO& eclIO() const; void writeInit(const std::function& map); void setTransmissibilities(const TransmissibilityType* globalTrans) { globalTrans_ = globalTrans; } void setSubStepReport(const SimulatorReportSingle& report) { sub_step_report_ = report; } void setSimulationReport(const SimulatorReport& report) { simulation_report_ = report; } protected: const TransmissibilityType& globalTrans() const; unsigned int gridEquilIdxToGridIdx(unsigned int elemIndex) const; void doWriteOutput(const int reportStepNum, const bool isSubStep, data::Solution&& localCellData, data::Wells&& localWellData, data::GroupAndNetworkValues&& localGroupAndNetworkData, data::Aquifers&& localAquiferData, WellTestState&& localWTestState, const Action::State& actionState, const UDQState& udqState, const SummaryState& summaryState, const std::vector& thresholdPressure, Scalar curTime, Scalar nextStepSize, bool doublePrecision); void evalSummary(int reportStepNum, Scalar curTime, const std::map& wbpData, const data::Wells& localWellData, const data::GroupAndNetworkValues& localGroupAndNetworkData, const std::map& localAquiferData, const std::map, double>& blockData, const std::map& miscSummaryData, const std::map>& regionData, const Inplace& inplace, const Inplace& initialInPlace, const EclInterRegFlowMap& interRegionFlowMap, SummaryState& summaryState, UDQState& udqState); CollectDataToIORankType collectToIORank_; const Grid& grid_; const GridView& gridView_; const Schedule& schedule_; const EclipseState& eclState_; const SummaryConfig& summaryConfig_; std::unique_ptr eclIO_; std::unique_ptr taskletRunner_; Scalar restartTimeStepSize_; const TransmissibilityType* globalTrans_ = nullptr; const Dune::CartesianIndexMapper& cartMapper_; const Dune::CartesianIndexMapper* equilCartMapper_; const EquilGrid* equilGrid_; std::vector wbp_index_list_; SimulatorReportSingle sub_step_report_; SimulatorReport simulation_report_; private: data::Solution computeTrans_(const std::unordered_map& cartesianToActive, const std::function& map) const; std::vector exportNncStructure_(const std::unordered_map& cartesianToActive, const std::function& map) const; }; } // namespace Opm #endif