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
31#include <opm/models/blackoil/blackoilmodel.hh>
32#include <opm/models/utils/propertysystem.hh>
33
34#include <opm/common/ErrorMacros.hpp>
35
37
38#include <stdexcept>
39
40namespace Opm
41{
42template <typename TypeTag>
43class FIBlackOilModel : public BlackOilModel<TypeTag>
44{
45 using ParentType = BlackOilModel<TypeTag>;
46 using Simulator = GetPropType<TypeTag, Properties::Simulator>;
47 using IntensiveQuantities = GetPropType<TypeTag, Properties::IntensiveQuantities>;
48 using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
49 using ThreadManager = GetPropType<TypeTag, Properties::ThreadManager>;
50 using GridView = GetPropType<TypeTag, Properties::GridView>;
51 using Element = typename GridView::template Codim<0>::Entity;
52 using ElementIterator = typename GridView::template Codim<0>::Iterator;
53 enum {
54 numEq = getPropValue<TypeTag, Properties::NumEq>(),
55 historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>(),
56 };
57
58public:
59 FIBlackOilModel(Simulator& simulator)
60 : BlackOilModel<TypeTag>(simulator)
61 {
62 }
63
64 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
65 {
66
67 this->invalidateIntensiveQuantitiesCache(timeIdx);
69 // loop over all elements...
70 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView_);
71#ifdef _OPENMP
72#pragma omp parallel
73#endif
74 {
75 ElementContext elemCtx(this->simulator_);
76 ElementIterator elemIt = threadedElemIt.beginParallel();
77 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
78 const Element& elem = *elemIt;
79 elemCtx.updatePrimaryStencil(elem);
80 elemCtx.updatePrimaryIntensiveQuantities(timeIdx);
81 }
82 }
83 OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantities: state error", this->simulator_.vanguard().grid().comm());
84 }
85
87 {
88 // loop over all elements
89 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView_);
91#ifdef _OPENMP
92#pragma omp parallel
93#endif
94 {
95 ElementContext elemCtx(this->simulator_);
96 auto elemIt = threadedElemIt.beginParallel();
97 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
98 if (elemIt->partitionType() != Dune::OverlapEntity) {
99 continue;
100 }
101 const Element& elem = *elemIt;
102 elemCtx.updatePrimaryStencil(elem);
103 // Mark cache for this element as invalid.
104 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
105 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
106 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
107 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
108 }
109 // Update for this element.
110 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
111 }
112 }
113 OPM_END_PARALLEL_TRY_CATCH("InvalideAndUpdateIntensiveQuantitiesOverlap: state error", this->simulator_.vanguard().grid().comm());
114 }
115
116 template <class GridSubDomain>
117 void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain& gridSubDomain) const
118 {
119 // loop over all elements in the subdomain
120 using GridViewType = decltype(gridSubDomain.view);
121 ThreadedEntityIterator<GridViewType, /*codim=*/0> threadedElemIt(gridSubDomain.view);
122#ifdef _OPENMP
123#pragma omp parallel
124#endif
125 {
126 ElementContext elemCtx(this->simulator_);
127 auto elemIt = threadedElemIt.beginParallel();
128 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
129 if (elemIt->partitionType() != Dune::InteriorEntity) {
130 continue;
131 }
132 const Element& elem = *elemIt;
133 elemCtx.updatePrimaryStencil(elem);
134 // Mark cache for this element as invalid.
135 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(timeIdx);
136 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
137 const unsigned globalIndex = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
138 this->setIntensiveQuantitiesCacheEntryValidity(globalIndex, timeIdx, false);
139 }
140 // Update for this element.
141 elemCtx.updatePrimaryIntensiveQuantities(/*timeIdx=*/0);
142 }
143 }
144 }
145
152 {
153 // Reset the current solution to the one of the
154 // previous time step so that we can start the next
155 // update at a physically meaningful solution.
156 // this->solution(/*timeIdx=*/0) = this->solution(/*timeIdx=*/1);
157 ParentType::updateFailed();
159 }
160
161 // standard flow
162 const IntensiveQuantities& intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
163 {
164 if (!this->enableIntensiveQuantityCache_) {
165 OPM_THROW(std::logic_error,
166 "Run without intensive quantites not enabled: Use --enable-intensive-quantity=true");
167 }
168 const auto* intquant = this->cachedIntensiveQuantities(globalIdx, timeIdx);
169 if (!intquant) {
170 OPM_THROW(std::logic_error, "Intensive quantites need to be updated in code");
171 }
172 return *intquant;
173 }
174};
175} // namespace Opm
176#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:172
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:138
Definition: FIBlackoilModel.hpp:44
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx, const GridSubDomain &gridSubDomain) const
Definition: FIBlackoilModel.hpp:117
void invalidateAndUpdateIntensiveQuantitiesOverlap(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:86
void updateFailed()
Called by the update() method if it was unsuccessful. This is primary a hook which the actual model c...
Definition: FIBlackoilModel.hpp:151
const IntensiveQuantities & intensiveQuantities(unsigned globalIdx, unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:162
void invalidateAndUpdateIntensiveQuantities(unsigned timeIdx) const
Definition: FIBlackoilModel.hpp:64
FIBlackOilModel(Simulator &simulator)
Definition: FIBlackoilModel.hpp:59
Definition: BlackoilPhases.hpp:27