pressureControl.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2017-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "pressureControl.H"
27 #include "findRefCell.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
32 (
33  const volScalarField& p,
34  const volScalarField& pRef,
35  const volScalarField& rho,
36  const dictionary& dict,
37  const bool pRefRequired
38 )
39 :
40  refCell_(-1),
41  refValue_(0),
42  pMax_("pMax", dimPressure, great),
43  pMin_("pMin", dimPressure, 0),
44  limitMaxP_(false),
45  limitMinP_(false)
46 {
47  bool pLimits = false;
48  scalar pMax = -great;
49  scalar pMin = great;
50 
51  // Set the reference cell and value for closed domain simulations
52  if (pRefRequired && setRefCell(p, pRef, dict, refCell_, refValue_))
53  {
54  pLimits = true;
55 
56  pMax = refValue_;
57  pMin = refValue_;
58  }
59 
60  if (dict.found("pMax") && dict.found("pMin"))
61  {
62  pMax_.value() = dict.lookup<scalar>("pMax");
63  limitMaxP_ = true;
64  pMin_.value() = dict.lookup<scalar>("pMin");
65  limitMinP_ = true;
66  }
67  else
68  {
69  const volScalarField::Boundary& pbf = p.boundaryField();
70  const volScalarField::Boundary& rhobf = rho.boundaryField();
71 
72  scalar rhoRefMax = -great;
73  scalar rhoRefMin = great;
74  bool rhoLimits = false;
75 
76  forAll(pbf, patchi)
77  {
78  if (pbf[patchi].fixesValue())
79  {
80  pLimits = true;
81  rhoLimits = true;
82 
83  pMax = max(pMax, max(pbf[patchi]));
84  pMin = min(pMin, min(pbf[patchi]));
85 
86  rhoRefMax = max(rhoRefMax, max(rhobf[patchi]));
87  rhoRefMin = min(rhoRefMin, min(rhobf[patchi]));
88  }
89  }
90 
91  reduce(rhoLimits, andOp<bool>());
92  if (rhoLimits)
93  {
94  reduce(pMax, maxOp<scalar>());
95  reduce(pMin, minOp<scalar>());
96 
97  reduce(rhoRefMax, maxOp<scalar>());
98  reduce(rhoRefMin, minOp<scalar>());
99  }
100 
101  if (dict.found("pMax"))
102  {
103  pMax_.value() = dict.lookup<scalar>("pMax");
104  limitMaxP_ = true;
105  }
106  else if (dict.found("pMaxFactor"))
107  {
108  if (!pLimits)
109  {
111  << "'pMaxFactor' specified rather than 'pMax'" << nl
112  << " but the corresponding reference pressure cannot"
113  " be evaluated from the boundary conditions." << nl
114  << " Please specify 'pMax' rather than 'pMaxFactor'"
115  << exit(FatalIOError);
116  }
117 
118  const scalar pMaxFactor(dict.lookup<scalar>("pMaxFactor"));
119  pMax_.value() = pMaxFactor*pMax;
120  limitMaxP_ = true;
121  }
122  else if (dict.found("rhoMax"))
123  {
124  // For backward-compatibility infer the pMax from rhoMax
125 
126  IOWarningInFunction(dict)
127  << "'rhoMax' specified rather than 'pMax' or 'pMaxFactor'"
128  << nl
129  << " This is supported for backward-compatibility but "
130  "'pMax' or 'pMaxFactor' are more reliable." << endl;
131 
132  if (!pLimits)
133  {
135  << "'rhoMax' specified rather than 'pMax'" << nl
136  << " but the corresponding reference pressure cannot"
137  " be evaluated from the boundary conditions." << nl
138  << " Please specify 'pMax' rather than 'rhoMax'"
139  << exit(FatalIOError);
140  }
141 
142  if (!rhoLimits)
143  {
145  << "'rhoMax' specified rather than 'pMaxFactor'" << nl
146  << " but the corresponding reference density cannot"
147  " be evaluated from the boundary conditions." << nl
148  << " Please specify 'pMaxFactor' rather than 'rhoMax'"
149  << exit(FatalIOError);
150  }
151 
152  dimensionedScalar rhoMax("rhoMax", dimDensity, dict);
153 
154  pMax_.value() = max(rhoMax.value()/rhoRefMax, 1)*pMax;
155  limitMaxP_ = true;
156  }
157 
158  if (dict.found("pMin"))
159  {
160  pMin_.value() = dict.lookup<scalar>("pMin");
161  limitMinP_ = true;
162  }
163  else if (dict.found("pMinFactor"))
164  {
165  if (!pLimits)
166  {
168  << "'pMinFactor' specified rather than 'pMin'" << nl
169  << " but the corresponding reference pressure cannot"
170  " be evaluated from the boundary conditions." << nl
171  << " Please specify 'pMin' rather than 'pMinFactor'"
172  << exit(FatalIOError);
173  }
174 
175  const scalar pMinFactor(dict.lookup<scalar>("pMinFactor"));
176  pMin_.value() = pMinFactor*pMin;
177  limitMinP_ = true;
178  }
179  else if (dict.found("rhoMin"))
180  {
181  // For backward-compatibility infer the pMin from rhoMin
182 
183  IOWarningInFunction(dict)
184  << "'rhoMin' specified rather than 'pMin' or 'pMinFactor'" << nl
185  << " This is supported for backward-compatibility but"
186  "'pMin' or 'pMinFactor' are more reliable." << endl;
187 
188  if (!pLimits)
189  {
191  << "'rhoMin' specified rather than 'pMin'" << nl
192  << " but the corresponding reference pressure cannot"
193  " be evaluated from the boundary conditions." << nl
194  << " Please specify 'pMin' rather than 'rhoMin'"
195  << exit(FatalIOError);
196  }
197 
198  if (!rhoLimits)
199  {
201  << "'rhoMin' specified rather than 'pMinFactor'" << nl
202  << " but the corresponding reference density cannot"
203  " be evaluated from the boundary conditions." << nl
204  << " Please specify 'pMinFactor' rather than 'rhoMin'"
205  << exit(FatalIOError);
206  }
207 
208  dimensionedScalar rhoMin("rhoMin", dimDensity, dict);
209 
210  pMin_.value() = min(rhoMin.value()/rhoRefMin, 1)*pMin;
211  limitMinP_ = true;
212  }
213  }
214 
215  if (limitMaxP_ || limitMinP_)
216  {
217  Info<< "pressureControl" << nl;
218 
219  if (limitMaxP_)
220  {
221  Info<< " pMax " << pMax_.value() << nl;
222  }
223 
224  if (limitMinP_)
225  {
226  Info<< " pMin " << pMin_.value() << nl;
227  }
228 
229  Info << endl;
230  }
231 }
232 
233 
235 (
236  const volScalarField& p,
237  const volScalarField& rho,
238  const dictionary& dict,
239  const bool pRefRequired
240 )
241 :
242  pressureControl(p, p, rho, dict, pRefRequired)
243 {}
244 
245 
246 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
247 
249 {
250  if (limitMaxP_ || limitMinP_)
251  {
252  if (limitMaxP_)
253  {
254  const scalar pMax = max(p).value();
255 
256  if (pMax > pMax_.value())
257  {
258  Info<< "pressureControl: p max " << pMax << endl;
259  p = min(p, pMax_);
260  }
261  }
262 
263  if (limitMinP_)
264  {
265  const scalar pMin = min(p).value();
266 
267  if (pMin < pMin_.value())
268  {
269  Info<< "pressureControl: p min " << pMin << endl;
270  p = max(p, pMin_);
271  }
272  }
273 
275 
276  return true;
277  }
278  else
279  {
280  return false;
281  }
282 }
283 
284 
285 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
pressureControl(const volScalarField &p, const volScalarField &pRef, const volScalarField &rho, const dictionary &dict, const bool pRefRequired=true)
Construct from the simple/pimple sub-dictionary.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Find the reference cell nearest (in index) to the given cell but which is not on a cyclic...
pressureControl & pressureControl
const Type & value() const
Return const reference to value.
bool limit(volScalarField &p) const
Limit the pressure if necessary and return true if so.
const dimensionSet dimPressure
dimensionedScalar pMin("pMin", dimPressure, fluid)
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimDensity
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
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:31
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
IOerror FatalIOError