FIBlackoilModel.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef FI_BLACK_OIL_MODEL_HPP
29#define FI_BLACK_OIL_MODEL_HPP
30
33
34#include <opm/common/ErrorMacros.hpp>
35
36#include <opm/grid/CpGrid.hpp>
37#include <opm/grid/utility/ElementChunks.hpp>
38
41
42#include <cstddef>
43#include <stdexcept>
44#include <type_traits>
45
46namespace Opm {
47
48template <typename TypeTag>
49class FIBlackOilModel : public BlackOilModel<TypeTag>
50{
57 using Element = typename GridView::template Codim<0>::Entity;
58 using ElementIterator = typename GridView::template Codim<0>::Iterator;
59 enum {
60 numEq = getPropValue<TypeTag, Properties::NumEq>(),
61 historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>(),
62 };
63
64 // The chunked and threaded iteration over elements in this class assumes that the number
65 // and order of elements is fixed, and is therefore constrained to only work with CpGrid.
66 // For example, ALUGrid supports refinement and can not assume that.
67 static constexpr bool gridIsUnchanging =
68 std::is_same_v<GetPropType<TypeTag, Properties::Grid>, Dune::CpGrid>;
69
70public:
71 explicit FIBlackOilModel(Simulator& simulator)
72 : BlackOilModel<TypeTag>(simulator)
73 , element_chunks_(this->gridView_,
74 Dune::Partitions::all,
75 ThreadManager::maxThreads())
76 {
77 }
78
79 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
80 {
81 this->invalidateIntensiveQuantitiesCache(timeIdx);
83 if constexpr (gridIsUnchanging) {
84#ifdef _OPENMP
85#pragma omp parallel for
86#endif
87 for (const auto& chunk : element_chunks_) {
88 ElementContext elemCtx(this->simulator_);
89 for (const auto& elem : chunk) {
90 elemCtx.updatePrimaryStencil(elem);
91 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
92 }
93 }
94 } else {
95 // Grid is possibly refined or otherwise changed between calls.
96 ElementContext elemCtx(this->simulator_);
97 for (const auto& elem : elements(this->gridView_)) {
98 elemCtx.updatePrimaryStencil(elem);
99 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
100 }
101 }
102 OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantities: state error",
103 this->simulator_.vanguard().grid().comm());
104 }
105
107 {
108 // loop over all elements
109 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView_);
111#ifdef _OPENMP
112#pragma omp parallel
113#endif
114 {
115 ElementContext elemCtx(this->simulator_);
116 auto elemIt = threadedElemIt.beginParallel();
117 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
118 if (elemIt->partitionType() != Dune::OverlapEntity) {
119 continue;
120 }
121 const Element& elem = *elemIt;
122 elemCtx.updatePrimaryStencil(elem);
123 // Mark cache for this element as invalid.
124 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
125 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
126 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
127 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
128 }
129 // Update for this element.
130 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
131 }
132 }
133 OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantitiesOverlap: state error",
134 this->simulator_.vanguard().grid().comm());
135 }
136
137 template <class GridSubDomain>
138 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain& gridSubDomain) const
139 {
140 // loop over all elements in the subdomain
141 using GridViewType = decltype(gridSubDomain.view);
142 ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(gridSubDomain.view);
143#ifdef _OPENMP
144#pragma omp parallel
145#endif
146 {
147 ElementContext elemCtx(this->simulator_);
148 auto elemIt = threadedElemIt.beginParallel();
149 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
150 if (elemIt->partitionType() != Dune::InteriorEntity) {
151 continue;
152 }
153 const Element& elem = *elemIt;
154 elemCtx.updatePrimaryStencil(elem);
155 // Mark cache for this element as invalid.
156 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
157 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
158 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
159 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
160 }
161 // Update for this element.
162 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
163 }
164 }
165 }
166
173 {
174 // Reset the current solution to the one of the
175 // previous time step so that we can start the next
176 // update at a physically meaningful solution.
177 // this->solution(/*timeIdx=*/0) = this->solution(/*timeIdx=*/1);
178 ParentType::updateFailed();
180 }
181
182 // standard flow
183 const IntensiveQuantities& intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
184 {
185 if (!this->enableIntensiveQuantityCache_) {
186 OPM_THROW(std::logic_error,
187 "Run without intensive quantites not enabled: "
188 "Use --enable-intensive-quantity=true");
189 }
190 const auto* intquant = this->cachedIntensiveQuantities(globalIdx, timeIdx);
191 if (!intquant) {
192 OPM_THROW(std::logic_error, "Intensive quantites need to be updated in code");
193 }
194 return *intquant;
195 }
196
197protected:
198 ElementChunks<GridView, Dune::Partitions::All> element_chunks_;
199};
200
201} // namespace Opm
202
203#endif // FI_BLACK_OIL_MODEL_HPP
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
A fully-implicit black-oil flow model.
Definition: blackoilmodel.hh:336
Definition: FIBlackoilModel.hpp:50
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain &gridSubDomain) const
Definition: FIBlackoilModel.hpp:138
void invalidateAndUpdateIntensiveQuantitiesOverlap(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:106
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: FIBlackoilModel.hpp:172
const IntensiveQuantities & intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:183
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:79
FIBlackOilModel(Simulator &simulator)
Definition: FIBlackoilModel.hpp:71
ElementChunks< GridView, Dune::Partitions::All > element_chunks_
Definition: FIBlackoilModel.hpp:198
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition: multiphasebasemodel.hh:168
Provides an STL-iterator like interface to iterate over the enties of a GridView in OpenMP threaded a...
Definition: threadedentityiterator.hh:42
bool isFinished(const EntityIterator &it) const
Definition: threadedentityiterator.hh:67
EntityIterator increment()
Definition: threadedentityiterator.hh:80
EntityIterator beginParallel()
Definition: threadedentityiterator.hh:54
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:39
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
The Opm property system, traits with inheritance.