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>
36
37#include <opm/grid/CpGrid.hpp>
38#include <opm/grid/utility/ElementChunks.hpp>
39
42
43#include <opm/material/fluidmatrixinteractions/EclMultiplexerMaterialParams.hpp>
44
45#include <cstddef>
46#include <stdexcept>
47#include <type_traits>
48
49namespace Opm {
50
51template <typename TypeTag>
52class FIBlackOilModel : public BlackOilModel<TypeTag>
53{
60 using Element = typename GridView::template Codim<0>::Entity;
61 using ElementIterator = typename GridView::template Codim<0>::Iterator;
62 enum {
63 numEq = getPropValue<TypeTag, Properties::NumEq>(),
64 historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>(),
65 };
66
67 // The chunked and threaded iteration over elements in this class assumes that the number
68 // and order of elements is fixed, and is therefore constrained to only work with CpGrid.
69 // For example, ALUGrid supports refinement and can not assume that.
70 static constexpr bool gridIsUnchanging =
71 std::is_same_v<GetPropType<TypeTag, Properties::Grid>, Dune::CpGrid>;
72
73 static constexpr bool avoidElementContext = getPropValue<TypeTag, Properties::AvoidElementContext>();
74
75public:
76 explicit FIBlackOilModel(Simulator& simulator)
77 : BlackOilModel<TypeTag>(simulator)
78 , element_chunks_(this->gridView_,
79 Dune::Partitions::all,
80 ThreadManager::maxThreads())
81 {
82 }
83
84 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
85 {
86 this->invalidateIntensiveQuantitiesCache(timeIdx);
87 if constexpr (gridIsUnchanging) {
88 if constexpr (avoidElementContext) {
89 updateCachedIntQuants(timeIdx);
90 return;
91 }
93#ifdef _OPENMP
94#pragma omp parallel for
95#endif
96 for (const auto& chunk : element_chunks_) {
97 ElementContext elemCtx(this->simulator_);
98 for (const auto& elem : chunk) {
99 elemCtx.updatePrimaryStencil(elem);
100 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
101 }
102 }
103 OPM_END_PARALLEL_TRY_CATCH("invalidateAndUpdateIntensiveQuantities: state error",
104 this->simulator_.vanguard().grid().comm());
105 } else {
106 // Grid is possibly refined or otherwise changed between calls.
107 ElementContext elemCtx(this->simulator_);
108 for (const auto& elem : elements(this->gridView_)) {
109 elemCtx.updatePrimaryStencil(elem);
110 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
111 }
112 }
113 }
114
116 {
117 // loop over all elements
118 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView_);
120#ifdef _OPENMP
121#pragma omp parallel
122#endif
123 {
124 ElementContext elemCtx(this->simulator_);
125 auto elemIt = threadedElemIt.beginParallel();
126 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
127 if (elemIt->partitionType() != Dune::OverlapEntity) {
128 continue;
129 }
130 const Element& elem = *elemIt;
131 elemCtx.updatePrimaryStencil(elem);
132 // Mark cache for this element as invalid.
133 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
134 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
135 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
136 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
137 }
138 // Update for this element.
139 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
140 }
141 }
142 OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantitiesOverlap: state error",
143 this->simulator_.vanguard().grid().comm());
144 }
145
146 template <class GridSubDomain>
147 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain& gridSubDomain) const
148 {
149 // loop over all elements in the subdomain
150 using GridViewType = decltype(gridSubDomain.view);
151 ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(gridSubDomain.view);
152#ifdef _OPENMP
153#pragma omp parallel
154#endif
155 {
156 ElementContext elemCtx(this->simulator_);
157 auto elemIt = threadedElemIt.beginParallel();
158 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
159 if (elemIt->partitionType() != Dune::InteriorEntity) {
160 continue;
161 }
162 const Element& elem = *elemIt;
163 elemCtx.updatePrimaryStencil(elem);
164 // Mark cache for this element as invalid.
165 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
166 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
167 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
168 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
169 }
170 // Update for this element.
171 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
172 }
173 }
174 }
175
182 {
183 // Reset the current solution to the one of the
184 // previous time step so that we can start the next
185 // update at a physically meaningful solution.
186 // this->solution(/*timeIdx=*/0) = this->solution(/*timeIdx=*/1);
187 ParentType::updateFailed();
189 }
190
191 // standard flow
192 const IntensiveQuantities& intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
193 {
194 if (!this->enableIntensiveQuantityCache_) {
195 OPM_THROW(std::logic_error,
196 "Run without intensive quantites not enabled: "
197 "Use --enable-intensive-quantity=true");
198 }
199
200 assert(timeIdx < this->cachedIntensiveQuantityHistorySize());
201 const auto* intquant = this->cachedIntensiveQuantities(globalIdx, timeIdx);
202 if (!intquant) {
203 OPM_THROW(std::logic_error, "Intensive quantites need to be updated in code");
204 }
205 return *intquant;
206 }
207
208protected:
209
210 template <EclMultiplexerApproach ApproachArg>
211 using EMD = EclMultiplexerDispatch<ApproachArg>;
212
213 void updateCachedIntQuants(const unsigned timeIdx) const
214 {
215 // Runtime dispatch on three-phase approach to compile-time fixed approach.
216 switch (this->simulator_.problem().materialLawManager()->threePhaseApproach()) {
217 case EclMultiplexerApproach::Stone1:
218 updateCachedIntQuants1<EclMultiplexerDispatch<EclMultiplexerApproach::Stone1>>(timeIdx);
219 break;
220
221 case EclMultiplexerApproach::Stone2:
222 updateCachedIntQuants1<EMD<EclMultiplexerApproach::Stone2>>(timeIdx);
223 break;
224
225 case EclMultiplexerApproach::Default:
226 updateCachedIntQuants1<EMD<EclMultiplexerApproach::Default>>(timeIdx);
227 break;
228
229 case EclMultiplexerApproach::TwoPhase:
230 updateCachedIntQuants1<EMD<EclMultiplexerApproach::TwoPhase>>(timeIdx);
231 break;
232
233 case EclMultiplexerApproach::OnePhase:
234 updateCachedIntQuants1<EMD<EclMultiplexerApproach::OnePhase>>(timeIdx);
235 break;
236 }
237 }
238
239 template <class EMDArg>
240 void updateCachedIntQuants1(const unsigned timeIdx) const
241 {
242 // Runtime dispatch on using piecewise linear saturation curves to compile-time fixed approach.
243 if (this->simulator_.problem().materialLawManager()->satCurveIsAllPiecewiseLinear()) {
244 using PL = SatCurveMultiplexerDispatch<SatCurveMultiplexerApproach::PiecewiseLinear>;
245 updateCachedIntQuantsLoop<EMDArg, PL>(timeIdx);
246 } else {
247 // TODO: Might want to set LET here, but need to check if partial use of LET is possible.
248 updateCachedIntQuantsLoop<EMDArg>(timeIdx);
249 }
250 }
251
252
253 template <class ...Args>
254 void updateCachedIntQuantsLoop(const unsigned timeIdx) const
255 {
256 const auto& elementMapper = this->simulator_.model().elementMapper();
257#ifdef _OPENMP
258#pragma omp parallel for
259#endif
260 for (const auto& chunk : element_chunks_) {
261 for (const auto& elem : chunk) {
262 this->template updateSingleCachedIntQuantUnchecked<Args...>(elementMapper.index(elem), timeIdx);
263 }
264 }
265 }
266
267 template <class ...Args>
268 void updateSingleCachedIntQuantUnchecked(const unsigned globalIdx, const unsigned timeIdx) const
269 {
270 // Get the cached data.
271 auto& intquant = this->intensiveQuantityCache_[timeIdx][globalIdx];
272 // Update it.
273 intquant.template update<Args...>(this->simulator_.problem(), this->solution(timeIdx)[globalIdx], globalIdx, timeIdx);
274 // Set the up-to-date flag.
275 this->intensiveQuantityCacheUpToDate_[timeIdx][globalIdx] = 1;
276 }
277
278 ElementChunks<GridView, Dune::Partitions::All> element_chunks_;
279};
280
281} // namespace Opm
282
283#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:53
void updateCachedIntQuantsLoop(const unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:254
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain &gridSubDomain) const
Definition: FIBlackoilModel.hpp:147
void updateCachedIntQuants(const unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:213
void updateSingleCachedIntQuantUnchecked(const unsigned globalIdx, const unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:268
void invalidateAndUpdateIntensiveQuantitiesOverlap(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:115
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: FIBlackoilModel.hpp:181
const IntensiveQuantities & intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:192
void updateCachedIntQuants1(const unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:240
EclMultiplexerDispatch< ApproachArg > EMD
Definition: FIBlackoilModel.hpp:211
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:84
FIBlackOilModel(Simulator &simulator)
Definition: FIBlackoilModel.hpp:76
ElementChunks< GridView, Dune::Partitions::All > element_chunks_
Definition: FIBlackoilModel.hpp:278
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.