// -*- 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. */ #ifndef EWOMS_ECL_GENERIC_WRITER_IMPL_HH #define EWOMS_ECL_GENERIC_WRITER_IMPL_HH #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #if HAVE_MPI #include #endif #if HAVE_MPI #include #endif #include #include #include #include #include #include #include #include #include #include namespace { /*! * \brief Detect whether two cells are direct vertical neighbours. * * I.e. have the same i and j index and all cartesian cells between them * along the vertical column are inactive. * * \tparam CM The type of the cartesian index mapper. * \param cartMapper The mapper onto cartesian indices. * \param cartesianToActive The mapping of cartesian indices to active indices. * \param smallGlobalIndex The cartesian cell index of the cell with smaller index * \param largeGlobalIndex The cartesian cell index of the cell with larger index * \return True if the cells have the same i and j indices and all cartesian cells * between them are inactive. */ bool directVerticalNeighbors(const std::array& cartDims, const std::unordered_map& cartesianToActive, int smallGlobalIndex, int largeGlobalIndex) { assert(smallGlobalIndex <= largeGlobalIndex); std::array ijk1, ijk2; auto globalToIjk = [cartDims](int gc) { std::array ijk; ijk[0] = gc % cartDims[0]; gc /= cartDims[0]; ijk[1] = gc % cartDims[1]; ijk[2] = gc / cartDims[1]; return ijk; }; ijk1 = globalToIjk(smallGlobalIndex); ijk2 = globalToIjk(largeGlobalIndex); assert(ijk2[2]>=ijk1[2]); if ( ijk1[0] == ijk2[0] && ijk1[1] == ijk2[1] && (ijk2[2] - ijk1[2]) > 1) { assert((largeGlobalIndex-smallGlobalIndex)%(cartDims[0]*cartDims[1])==0); for ( int gi = smallGlobalIndex + cartDims[0] * cartDims[1]; gi < largeGlobalIndex; gi += cartDims[0] * cartDims[1] ) { if ( cartesianToActive.find( gi ) != cartesianToActive.end() ) { return false; } } return true; } else return false; } std::unordered_map getInterRegFlowsAsMap(const Opm::EclInterRegFlowMap& map) { auto maps = std::unordered_map{}; const auto& regionNames = map.names(); auto flows = map.getInterRegFlows(); const auto nmap = regionNames.size(); maps.reserve(nmap); for (auto mapID = 0*nmap; mapID < nmap; ++mapID) { maps.emplace(regionNames[mapID], std::move(flows[mapID])); } return maps; } struct EclWriteTasklet : public Opm::TaskletInterface { Opm::Action::State actionState_; Opm::WellTestState wtestState_; Opm::SummaryState summaryState_; Opm::UDQState udqState_; Opm::EclipseIO& eclIO_; int reportStepNum_; bool isSubStep_; double secondsElapsed_; Opm::RestartValue restartValue_; bool writeDoublePrecision_; explicit EclWriteTasklet(const Opm::Action::State& actionState, const Opm::WellTestState& wtestState, const Opm::SummaryState& summaryState, const Opm::UDQState& udqState, Opm::EclipseIO& eclIO, int reportStepNum, bool isSubStep, double secondsElapsed, Opm::RestartValue restartValue, bool writeDoublePrecision) : actionState_(actionState) , wtestState_(wtestState) , summaryState_(summaryState) , udqState_(udqState) , eclIO_(eclIO) , reportStepNum_(reportStepNum) , isSubStep_(isSubStep) , secondsElapsed_(secondsElapsed) , restartValue_(std::move(restartValue)) , writeDoublePrecision_(writeDoublePrecision) {} // callback to eclIO serial writeTimeStep method void run() { this->eclIO_.writeTimeStep(this->actionState_, this->wtestState_, this->summaryState_, this->udqState_, this->reportStepNum_, this->isSubStep_, this->secondsElapsed_, std::move(this->restartValue_), this->writeDoublePrecision_); } }; } namespace Opm { template EclGenericWriter:: 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 ) : collectToIORank_(grid, equilGrid, gridView, cartMapper, equilCartMapper, summaryConfig.fip_regions_interreg_flow()) , grid_ (grid) , gridView_ (gridView) , schedule_ (schedule) , eclState_ (eclState) , summaryConfig_ (summaryConfig) , cartMapper_ (cartMapper) , equilCartMapper_(equilCartMapper) , equilGrid_ (equilGrid) { if (this->collectToIORank_.isIORank()) { this->eclIO_ = std::make_unique (this->eclState_, UgGridHelpers::createEclipseGrid(*equilGrid, eclState_.getInputGrid()), this->schedule_, this->summaryConfig_, "", enableEsmry); } // create output thread if enabled and rank is I/O rank // async output is enabled by default if pthread are enabled int numWorkerThreads = 0; if (enableAsyncOutput && collectToIORank_.isIORank()) { numWorkerThreads = 1; } this->taskletRunner_.reset(new TaskletRunner(numWorkerThreads)); } template const EclipseIO& EclGenericWriter:: eclIO() const { assert(eclIO_); return *eclIO_; } template void EclGenericWriter:: writeInit(const std::function& map) { if (collectToIORank_.isIORank()) { std::map> integerVectors; if (collectToIORank_.isParallel()) { integerVectors.emplace("MPI_RANK", collectToIORank_.globalRanks()); } auto cartMap = cartesianToCompressed(equilGrid_->size(0), UgGridHelpers::globalCell(*equilGrid_)); eclIO_->writeInitial(computeTrans_(cartMap, map), integerVectors, exportNncStructure_(cartMap, map)); } #if HAVE_MPI if (collectToIORank_.isParallel()) { const auto& comm = grid_.comm(); Opm::EclMpiSerializer ser(comm); ser.broadcast(outputNnc_); } #endif } template data::Solution EclGenericWriter:: computeTrans_(const std::unordered_map& cartesianToActive, const std::function& map) const { const auto& cartMapper = *equilCartMapper_; const auto& cartDims = cartMapper.cartesianDimensions(); const int globalSize = cartDims[0] * cartDims[1] * cartDims[2]; auto tranx = data::CellData { UnitSystem::measure::transmissibility, std::vector(globalSize, 0.0), data::TargetType::INIT }; auto trany = tranx; auto tranz = tranx; using GlobalGridView = typename EquilGrid::LeafGridView; const GlobalGridView& globalGridView = equilGrid_->leafGridView(); using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper; GlobElementMapper globalElemMapper(globalGridView, Dune::mcmgElementLayout()); for (const auto& elem : elements(globalGridView)) { for (const auto& is : intersections(globalGridView, elem)) { if (!is.neighbor()) continue; // intersection is on the domain boundary unsigned c1 = globalElemMapper.index(is.inside()); unsigned c2 = globalElemMapper.index(is.outside()); if (c1 > c2) continue; // we only need to handle each connection once, thank you. // Ordering of compressed and uncompressed index should be the same const int cartIdx1 = cartMapper.cartesianIndex( c1 ); const int cartIdx2 = cartMapper.cartesianIndex( c2 ); // Ordering of compressed and uncompressed index should be the same assert(cartIdx1 <= cartIdx2); int gc1 = std::min(cartIdx1, cartIdx2); int gc2 = std::max(cartIdx1, cartIdx2); // Re-ordering in case of non-empty mapping between equilGrid to grid if (map) { c1 = map(c1); // equilGridToGrid map c2 = map(c2); } if (gc2 - gc1 == 1 && cartDims[0] > 1 ) { tranx.data[gc1] = globalTrans().transmissibility(c1, c2); continue; // skip other if clauses as they are false, last one needs some computation } if (gc2 - gc1 == cartDims[0] && cartDims[1] > 1) { trany.data[gc1] = globalTrans().transmissibility(c1, c2); continue; // skipt next if clause as it needs some computation } if ( gc2 - gc1 == cartDims[0]*cartDims[1] || directVerticalNeighbors(cartDims, cartesianToActive, gc1, gc2)) tranz.data[gc1] = globalTrans().transmissibility(c1, c2); } } return { {"TRANX", tranx}, {"TRANY", trany}, {"TRANZ", tranz}, }; } template std::vector EclGenericWriter:: exportNncStructure_(const std::unordered_map& cartesianToActive, const std::function& map) const { std::size_t nx = eclState_.getInputGrid().getNX(); std::size_t ny = eclState_.getInputGrid().getNY(); auto nncData = eclState_.getInputNNC().input(); const auto& unitSystem = eclState_.getDeckUnitSystem(); for( const auto& entry : nncData ) { // test whether NNC is not a neighboring connection // cell2>=cell1 holds due to sortNncAndApplyEditnnc assert( entry.cell2 >= entry.cell1 ); auto cellDiff = entry.cell2 - entry.cell1; if (cellDiff != 1 && cellDiff != nx && cellDiff != nx*ny) { auto tt = unitSystem.from_si(UnitSystem::measure::transmissibility, entry.trans); // Eclipse ignores NNCs (with EDITNNC applied) that are small. Seems like the threshold is 1.0e-6 if ( tt >= 1.0e-6 ) outputNnc_.emplace_back(entry.cell1, entry.cell2, entry.trans); } } using GlobalGridView = typename EquilGrid::LeafGridView; const GlobalGridView& globalGridView = equilGrid_->leafGridView(); using GlobElementMapper = Dune::MultipleCodimMultipleGeomTypeMapper; GlobElementMapper globalElemMapper(globalGridView, Dune::mcmgElementLayout()); // Cartesian index mapper for the serial I/O grid const auto& equilCartMapper = *equilCartMapper_; const auto& cartDims = cartMapper_.cartesianDimensions(); for (const auto& elem : elements(globalGridView)) { for (const auto& is : intersections(globalGridView, elem)) { if (!is.neighbor()) continue; // intersection is on the domain boundary unsigned c1 = globalElemMapper.index(is.inside()); unsigned c2 = globalElemMapper.index(is.outside()); if (c1 > c2) continue; // we only need to handle each connection once, thank you. std::size_t cc1 = equilCartMapper.cartesianIndex( c1 ); std::size_t cc2 = equilCartMapper.cartesianIndex( c2 ); if ( cc2 < cc1 ) std::swap(cc1, cc2); auto cellDiff = cc2 - cc1; // Re-ordering in case of non-empty mapping between equilGrid to grid if (map) { c1 = map(c1); // equilGridToGrid map c2 = map(c2); } if (cellDiff != 1 && cellDiff != nx && cellDiff != nx*ny && !directVerticalNeighbors(cartDims, cartesianToActive, cc1, cc2)) { // We need to check whether an NNC for this face was also specified // via the NNC keyword in the deck (i.e. in the first origNncSize entries. auto t = globalTrans().transmissibility(c1, c2); auto candidate = std::lower_bound(nncData.begin(), nncData.end(), NNCdata(cc1, cc2, 0.0)); while ( candidate != nncData.end() && candidate->cell1 == cc1 && candidate->cell2 == cc2) { t -= candidate->trans; ++candidate; } // eclipse ignores NNCs with zero transmissibility (different threshold than for NNC // with corresponding EDITNNC above). In addition we do set small transmissibilties // to zero when setting up the simulator. These will be ignored here, too. auto tt = unitSystem.from_si(UnitSystem::measure::transmissibility, std::abs(t)); if ( tt > 1e-12 ) outputNnc_.push_back({cc1, cc2, t}); } } } return outputNnc_; } template void EclGenericWriter:: 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, bool isFlowsn, std::array, std::vector>>, 3>&& flowsn, bool isFloresn, std::array, std::vector>>, 3>&& floresn) { const auto isParallel = this->collectToIORank_.isParallel(); const bool needsReordering = this->collectToIORank_.doesNeedReordering(); RestartValue restartValue { (isParallel || needsReordering) ? this->collectToIORank_.globalCellData() : std::move(localCellData), isParallel ? this->collectToIORank_.globalWellData() : std::move(localWellData), isParallel ? this->collectToIORank_.globalGroupAndNetworkData() : std::move(localGroupAndNetworkData), isParallel ? this->collectToIORank_.globalAquiferData() : std::move(localAquiferData) }; if (eclState_.getSimulationConfig().useThresholdPressure()) { restartValue.addExtra("THRESHPR", UnitSystem::measure::pressure, thresholdPressure); } // Add suggested next timestep to extra data. if (! isSubStep) { restartValue.addExtra("OPMEXTRA", std::vector(1, nextStepSize)); } // Add nnc flows and flores. if (isFlowsn) { const auto flowsn_global = isParallel ? this->collectToIORank_.globalFlowsn() : std::move(flowsn); for (const auto& flows : flowsn_global) { if (flows.first.empty()) continue; if (flows.first == "FLOGASN+") { restartValue.addExtra(flows.first, UnitSystem::measure::gas_surface_rate, flows.second.second); } else { restartValue.addExtra(flows.first, UnitSystem::measure::liquid_surface_rate, flows.second.second); } } } if (isFloresn) { const auto floresn_global = isParallel ? this->collectToIORank_.globalFloresn() : std::move(floresn); for (const auto& flores : floresn_global) { if (flores.first.empty()) continue; restartValue.addExtra(flores.first, UnitSystem::measure::rate, flores.second.second); } } // first, create a tasklet to write the data for the current time // step to disk auto eclWriteTasklet = std::make_shared( actionState, isParallel ? this->collectToIORank_.globalWellTestState() : std::move(localWTestState), summaryState, udqState, *this->eclIO_, reportStepNum, isSubStep, curTime, std::move(restartValue), doublePrecision); // then, make sure that the previous I/O request has been completed // and the number of incomplete tasklets does not increase between // time steps this->taskletRunner_->barrier(); // finally, start a new output writing job this->taskletRunner_->dispatch(std::move(eclWriteTasklet)); } template void EclGenericWriter:: evalSummary(const int reportStepNum, const Scalar curTime, const data::Wells& localWellData, const data::WellBlockAveragePressures& localWBPData, 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& interRegFlows, SummaryState& summaryState, UDQState& udqState) { if (collectToIORank_.isIORank()) { const auto& summary = eclIO_->summary(); const auto& wellData = this->collectToIORank_.isParallel() ? this->collectToIORank_.globalWellData() : localWellData; const auto& wbpData = this->collectToIORank_.isParallel() ? this->collectToIORank_.globalWBPData() : localWBPData; const auto& groupAndNetworkData = this->collectToIORank_.isParallel() ? this->collectToIORank_.globalGroupAndNetworkData() : localGroupAndNetworkData; const auto& aquiferData = this->collectToIORank_.isParallel() ? this->collectToIORank_.globalAquiferData() : localAquiferData; summary.eval(summaryState, reportStepNum, curTime, wellData, wbpData, groupAndNetworkData, miscSummaryData, initialInPlace, inplace, regionData, blockData, aquiferData, getInterRegFlowsAsMap(interRegFlows)); // Off-by-one-fun: The reportStepNum argument corresponds to the // report step these results will be written to, whereas the // argument to UDQ function evaluation corresponds to the report // step we are currently on. const auto udq_step = reportStepNum - 1; this->schedule_.getUDQConfig(udq_step) .eval(udq_step, this->schedule_, this->schedule_.wellMatcher(udq_step), summaryState, udqState); } #if HAVE_MPI if (collectToIORank_.isParallel()) { EclMpiSerializer ser(grid_.comm()); ser.append(summaryState); } #endif } template const typename EclGenericWriter::TransmissibilityType& EclGenericWriter:: globalTrans() const { assert (globalTrans_); return *globalTrans_; } } // namespace Opm #endif