pressureControl.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "pressureControl.H"
29#include "findRefCell.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
34(
35 const volScalarField& p,
36 const volScalarField& rho,
37 const dictionary& dict,
38 const bool pRefRequired
39)
40:
41 refCell_(-1),
42 refValue_(0),
43 pMax_("pMax", dimPressure, GREAT),
44 pMin_("pMin", dimPressure, Zero),
45 limitMaxP_(false),
46 limitMinP_(false)
47{
48 bool pLimits = false;
49 scalar pMax = -GREAT;
50 scalar pMin = GREAT;
51
52 // Set the reference cell and value for closed domain simulations
53 if (pRefRequired && setRefCell(p, dict, refCell_, refValue_))
54 {
55 pLimits = true;
56
57 pMax = refValue_;
58 pMin = refValue_;
59 }
60
61 if (dict.found("pMax") && dict.found("pMin"))
62 {
63 dict.readEntry("pMax", pMax_.value()); limitMaxP_ = true;
64 dict.readEntry("pMin", pMin_.value()); limitMinP_ = true;
65 }
66 else
67 {
69 const volScalarField::Boundary& rhobf = rho.boundaryField();
70
71 scalar rhoRefMax = -GREAT;
72 scalar rhoRefMin = GREAT;
73 bool rhoLimits = false;
74
75 forAll(pbf, patchi)
76 {
77 if (pbf[patchi].fixesValue())
78 {
79 pLimits = true;
80 rhoLimits = true;
81
82 pMax = max(pMax, max(pbf[patchi]));
83 pMin = min(pMin, min(pbf[patchi]));
84
85 rhoRefMax = max(rhoRefMax, max(rhobf[patchi]));
86 rhoRefMin = min(rhoRefMin, min(rhobf[patchi]));
87 }
88 }
89
90 reduce(rhoLimits, andOp<bool>());
91 if (rhoLimits)
92 {
93 reduce(pMax, maxOp<scalar>());
95
96 reduce(rhoRefMax, maxOp<scalar>());
97 reduce(rhoRefMin, minOp<scalar>());
98 }
99
100 if (dict.readIfPresent("pMax", pMax_.value()))
101 {
102 limitMaxP_ = true;
103 }
104 else if (dict.found("pMaxFactor"))
105 {
106 if (!pLimits)
107 {
109 << "'pMaxFactor' specified rather than 'pMax'" << nl
110 << " but the corresponding reference pressure cannot"
111 " be evaluated from the boundary conditions." << nl
112 << " Please specify 'pMax' rather than 'pMaxFactor'"
113 << exit(FatalIOError);
114 }
115
116 pMax_.value() = pMax * dict.get<scalar>("pMaxFactor");
117 limitMaxP_ = true;
118 }
119 else if (dict.found("rhoMax"))
120 {
121 // For backward-compatibility infer the pMax from rhoMax
122
124 << "'rhoMax' specified rather than 'pMax' or 'pMaxFactor'"
125 << nl
126 << " This is supported for backward-compatibility but"
127 " 'pMax' or 'pMaxFactor' are more reliable." << endl;
128
129 if (!pLimits)
130 {
132 << "'rhoMax' specified rather than 'pMax'" << nl
133 << " but the corresponding reference pressure cannot"
134 " be evaluated from the boundary conditions." << nl
135 << " Please specify 'pMax' rather than 'rhoMax'"
136 << exit(FatalIOError);
137 }
138
139 if (!rhoLimits)
140 {
142 << "'rhoMax' specified rather than 'pMaxFactor'" << nl
143 << " but the corresponding reference density cannot"
144 " be evaluated from the boundary conditions." << nl
145 << " Please specify 'pMaxFactor' rather than 'rhoMax'"
146 << exit(FatalIOError);
147 }
148
150
151 pMax_.value() = max(rhoMax.value()/rhoRefMax, 1)*pMax;
152 limitMaxP_ = true;
153 }
154
155 if (dict.readIfPresent("pMin", pMin_.value()))
156 {
157 limitMinP_ = true;
158 }
159 else if (dict.found("pMinFactor"))
160 {
161 if (!pLimits)
162 {
164 << "'pMinFactor' specified rather than 'pMin'" << nl
165 << " but the corresponding reference pressure cannot"
166 " be evaluated from the boundary conditions." << nl
167 << " Please specify 'pMin' rather than 'pMinFactor'"
168 << exit(FatalIOError);
169 }
170
171 pMin_.value() = pMin * dict.get<scalar>("pMinFactor");
172 limitMinP_ = true;
173 }
174 else if (dict.found("rhoMin"))
175 {
176 // For backward-compatibility infer the pMin from rhoMin
177
179 << "'rhoMin' specified rather than 'pMin' or 'pMinFactor'" << nl
180 << " This is supported for backward-compatibility but"
181 " 'pMin' or 'pMinFactor' are more reliable." << endl;
182
183 if (!pLimits)
184 {
186 << "'rhoMin' specified rather than 'pMin'" << nl
187 << " but the corresponding reference pressure cannot"
188 " be evaluated from the boundary conditions." << nl
189 << " Please specify 'pMin' rather than 'rhoMin'"
190 << exit(FatalIOError);
191 }
192
193 if (!rhoLimits)
194 {
196 << "'rhoMin' specified rather than 'pMinFactor'" << nl
197 << " but the corresponding reference density cannot"
198 " be evaluated from the boundary conditions." << nl
199 << " Please specify 'pMinFactor' rather than 'rhoMin'"
200 << exit(FatalIOError);
201 }
202
204
205 pMin_.value() = min(rhoMin.value()/rhoRefMin, 1)*pMin;
206 limitMinP_ = true;
207 }
208 }
209
210 if (limitMaxP_ || limitMinP_)
211 {
212 Info<< "pressureControl" << nl;
213
214 if (limitMaxP_)
215 {
216 Info<< " pMax " << pMax_.value() << nl;
217 }
218
219 if (limitMinP_)
220 {
221 Info<< " pMin " << pMin_.value() << nl;
222 }
223
224 Info << endl;
225 }
226}
227
228
229// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
230
232{
233 if (limitMaxP_ || limitMinP_)
234 {
235 if (limitMaxP_)
236 {
237 const scalar pMax = max(p).value();
238
239 if (pMax > pMax_.value())
240 {
241 Info<< "pressureControl: p max " << pMax << endl;
242 p = min(p, pMax_);
243 }
244 }
245
246 if (limitMinP_)
247 {
248 const scalar pMin = min(p).value();
249
250 if (pMin < pMin_.value())
251 {
252 Info<< "pressureControl: p min " << pMin << endl;
253 p = max(p, pMin_);
254 }
255 }
256
257 return true;
258 }
259
260 return false;
261}
262
263
264// ************************************************************************* //
const dimensionedScalar & pMin
const Boundary & boundaryField() const
Return const-reference to the boundary field.
friend complex limit(const complex &c1, const complex &c2)
Definition: complexI.H:263
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
const Type & value() const
Return const reference to value.
Provides controls for the pressure reference is closed-volume simulations and a general method for li...
volScalarField & p
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic,...
const dimensionedScalar rhoMin
const dimensionedScalar rhoMax
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const dimensionSet dimPressure
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
bool setRefCell(const volScalarField &field, const volScalarField &fieldRef, const dictionary &dict, label &refCelli, scalar &refValue, const bool forceReference=false)
If the field fieldRef needs referencing find the reference cell nearest.
Definition: findRefCell.C:34
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333