// -*- 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::FIBlackOilModel */ #ifndef FI_BLACK_OIL_MODEL_HPP #define FI_BLACK_OIL_MODEL_HPP #include #include #include #include #include namespace Opm { template class FIBlackOilModel : public BlackOilModel { using ParentType = BlackOilModel; using Simulator = GetPropType; using IntensiveQuantities = GetPropType; using ElementContext = GetPropType; using ThreadManager = GetPropType; using GridView = GetPropType; using Element = typename GridView::template Codim<0>::Entity; using ElementIterator = typename GridView::template Codim<0>::Iterator; enum { numEq = getPropValue(), historySize = getPropValue(), }; public: FIBlackOilModel(Simulator& simulator) : BlackOilModel(simulator) { } void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const { this->invalidateIntensiveQuantitiesCache(timeIdx); OPM_BEGIN_PARALLEL_TRY_CATCH() // loop over all elements... ThreadedEntityIterator threadedElemIt(this->gridView_); #ifdef _OPENMP #pragma omp parallel #endif { ElementContext elemCtx(this->simulator_); ElementIterator elemIt = threadedElemIt.beginParallel(); for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) { const Element& elem = *elemIt; elemCtx.updatePrimaryStencil(elem); elemCtx.updatePrimaryIntensiveQuantities(timeIdx); } } OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantities: state error", this->simulator_.vanguard().grid().comm()); } void invalidateAndUpdateIntensiveQuantitiesOverlap(unsigned timeIdx) const { // loop over all elements ThreadedEntityIterator threadedElemIt(this->gridView_); OPM_BEGIN_PARALLEL_TRY_CATCH() #ifdef _OPENMP #pragma omp parallel #endif { ElementContext elemCtx(this->simulator_); auto elemIt = threadedElemIt.beginParallel(); for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) { if (elemIt->partitionType() != Dune::OverlapEntity) { continue; } const Element& elem = *elemIt; elemCtx.updatePrimaryStencil(elem); // Mark cache for this element as invalid. const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx); for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) { const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx); this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false); } // Update for this element. elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0); } } OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantitiesOverlap: state error", this->simulator_.vanguard().grid().comm()); } template void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain& gridSubDomain) const { // loop over all elements in the subdomain using GridViewType = decltype(gridSubDomain.view); ThreadedEntityIterator threadedElemIt(gridSubDomain.view); #ifdef _OPENMP #pragma omp parallel #endif { ElementContext elemCtx(this->simulator_); auto elemIt = threadedElemIt.beginParallel(); for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) { if (elemIt->partitionType() != Dune::InteriorEntity) { continue; } const Element& elem = *elemIt; elemCtx.updatePrimaryStencil(elem); // Mark cache for this element as invalid. const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx); for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) { const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx); this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false); } // Update for this element. elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0); } } } /*! * \brief Called by the update() method if it was * unsuccessful. This is primary a hook which the actual * model can overload. */ void updateFailed() { // Reset the current solution to the one of the // previous time step so that we can start the next // update at a physically meaningful solution. // this->solution(/*timeIdx=*/0) = this->solution(/*timeIdx=*/1); ParentType::updateFailed(); invalidateAndUpdateIntensiveQuantities(/*timeIdx=*/0); } // standard flow const IntensiveQuantities& intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const { if (!this->enableIntensiveQuantityCache_) { OPM_THROW(std::logic_error, "Run without intensive quantites not enabled: Use --enable-intensive-quantity=true"); } const auto* intquant = this->cachedIntensiveQuantities(globalIdx, timeIdx); if (!intquant) { OPM_THROW(std::logic_error, "Intensive quantites need to be updated in code"); } return *intquant; } }; } // namespace Opm #endif // FI_BLACK_OIL_MODEL_HPP