FlowThresholdPressure.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 OPM_FLOW_THRESHOLD_PRESSURE_HPP
29#define OPM_FLOW_THRESHOLD_PRESSURE_HPP
30
31#include <opm/material/densead/Evaluation.hpp>
32#include <opm/material/densead/Math.hpp>
33
37
39
40#include <algorithm>
41
42namespace Opm {
43
54template <class TypeTag>
55class FlowThresholdPressure : public GenericThresholdPressure<GetPropType<TypeTag, Properties::Grid>,
56 GetPropType<TypeTag, Properties::GridView>,
57 GetPropType<TypeTag, Properties::ElementMapper>,
58 GetPropType<TypeTag, Properties::Scalar>>
59{
69
70 enum { numPhases = FluidSystem::numPhases };
71
72public:
73 explicit FlowThresholdPressure(const Simulator& simulator)
74 : BaseType(simulator.vanguard().cartesianIndexMapper(),
75 simulator.vanguard().gridView(),
76 simulator.model().elementMapper(),
77 simulator.vanguard().eclState())
78 , simulator_(simulator)
79 {
80 }
81
86 {
88 if (this->enableThresholdPressure_ && !this->thpresDefault_.empty() && !this->restart_) {
89 this->computeDefaultThresholdPressures_();
91 simulator_.vanguard().gridView().comm().max(&this->thpres_[0], this->thpres_.size());
92 }
93 }
94
95private:
96 // compute the defaults of the threshold pressures using the initial condition
97 void computeDefaultThresholdPressures_()
98 {
99 const auto& vanguard = simulator_.vanguard();
100 const auto& gridView = vanguard.gridView();
101
102 using Toolbox = MathToolbox<Evaluation>;
103 // loop over the whole grid and compute the maximum gravity adjusted pressure
104 // difference between two EQUIL regions.
105 ElementContext elemCtx(simulator_);
106 for (const auto& elem : elements(gridView, Dune::Partitions::interior)) {
107 elemCtx.updateAll(elem);
108 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
109
110 for (unsigned scvfIdx = 0; scvfIdx < stencil.numInteriorFaces(); ++ scvfIdx) {
111 const auto& face = stencil.interiorFace(scvfIdx);
112
113 unsigned i = face.interiorIndex();
114 unsigned j = face.exteriorIndex();
115
116 unsigned insideElemIdx = elemCtx.globalSpaceIndex(i, /*timeIdx=*/0);
117 unsigned outsideElemIdx = elemCtx.globalSpaceIndex(j, /*timeIdx=*/0);
118
119 unsigned equilRegionInside = this->elemEquilRegion_[insideElemIdx];
120 unsigned equilRegionOutside = this->elemEquilRegion_[outsideElemIdx];
121
122 if (equilRegionInside == equilRegionOutside)
123 // the current face is not at the boundary between EQUIL regions!
124 continue;
125
126 // don't include connections with negligible flow
127 const Evaluation& trans = simulator_.problem().transmissibility(elemCtx, i, j);
128 Scalar faceArea = face.area();
129 if (std::abs(faceArea*getValue(trans)) < 1e-18)
130 continue;
131
132 // determine the maximum difference of the pressure of any phase over the
133 // intersection
134 Scalar pth = 0.0;
135 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, /*timeIdx=*/0);
136 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
137 unsigned upIdx = extQuants.upstreamIndex(phaseIdx);
138 const auto& up = elemCtx.intensiveQuantities(upIdx, /*timeIdx=*/0);
139
140 if (up.mobility(phaseIdx) > 0.0) {
141 Scalar phaseVal = Toolbox::value(extQuants.pressureDifference(phaseIdx));
142 pth = std::max(pth, std::abs(phaseVal));
143 }
144 }
145
146 int offset1 = equilRegionInside*this->numEquilRegions_ + equilRegionOutside;
147 int offset2 = equilRegionOutside*this->numEquilRegions_ + equilRegionInside;
148
149 this->thpresDefault_[offset1] = std::max(this->thpresDefault_[offset1], pth);
150 this->thpresDefault_[offset2] = std::max(this->thpresDefault_[offset2], pth);
151 }
152 }
153
154 // make sure that the threshold pressures is consistent for parallel
155 // runs. (i.e. take the maximum of all processes)
156 for (unsigned i = 0; i < this->thpresDefault_.size(); ++i)
157 this->thpresDefault_[i] = gridView.comm().max(this->thpresDefault_[i]);
158
159 this->logPressures();
160 }
161
162 const Simulator& simulator_;
163};
164
165} // namespace Opm
166
167#endif // OPM_FLOW_THRESHOLD_PRESSURE_HPP
This class calculates the threshold pressure for grid faces according to the Eclipse Reference Manual...
Definition: FlowThresholdPressure.hpp:59
FlowThresholdPressure(const Simulator &simulator)
Definition: FlowThresholdPressure.hpp:73
void finishInit()
Actually compute the threshold pressures over a face as a pre-compute step.
Definition: FlowThresholdPressure.hpp:85
Definition: GenericThresholdPressure.hpp:43
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common properties required by the porous medium multi-phase models.
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.