limitPressure.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) 2021 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 "limitPressure.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace fv
35 {
36  defineTypeNameAndDebug(limitPressure, 0);
38  (
39  fvConstraint,
40  limitPressure,
41  dictionary
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::limitPressure::readCoeffs()
50 {
51  const dictionary& dict(coeffs());
52 
53  if (dict.found("min") && dict.found("max"))
54  {
55  pMin_.value() = dict.lookup<scalar>("min");
56  limitMinP_ = true;
57 
58  pMax_.value() = dict.lookup<scalar>("max");
59  limitMaxP_ = true;
60  }
61  else
62  {
63  const volScalarField& p = mesh().lookupObject<volScalarField>("p");
64  const volScalarField::Boundary& pbf = p.boundaryField();
65 
66  bool pLimits = false;
67  scalar pMin = vGreat;
68  scalar pMax = -vGreat;
69 
70  forAll(pbf, patchi)
71  {
72  if (pbf[patchi].fixesValue())
73  {
74  pLimits = true;
75 
76  pMin = min(pMin, min(pbf[patchi]));
77  pMax = max(pMax, max(pbf[patchi]));
78  }
79  }
80 
81  reduce(pLimits, andOp<bool>());
82  if (pLimits)
83  {
84  reduce(pMax, maxOp<scalar>());
85  reduce(pMin, minOp<scalar>());
86  }
87 
88  if (dict.found("min"))
89  {
90  pMin_.value() = dict.lookup<scalar>("min");
91  limitMinP_ = true;
92  }
93  else if (dict.found("minFactor"))
94  {
95  if (!pLimits)
96  {
98  << "'minFactor' specified rather than 'min'" << nl
99  << " but the corresponding reference pressure cannot"
100  " be evaluated from the boundary conditions." << nl
101  << " Please specify 'min' rather than 'minFactor'"
102  << exit(FatalIOError);
103  }
104 
105  const scalar pMinFactor(dict.lookup<scalar>("minFactor"));
106  pMin_.value() = pMinFactor*pMin;
107  limitMinP_ = true;
108  }
109 
110  if (dict.found("max"))
111  {
112  pMax_.value() = dict.lookup<scalar>("max");
113  limitMaxP_ = true;
114  }
115  else if (dict.found("maxFactor"))
116  {
117  if (!pLimits)
118  {
120  << "'maxFactor' specified rather than 'max'" << nl
121  << " but the corresponding reference pressure cannot"
122  " be evaluated from the boundary conditions." << nl
123  << " Please specify 'max' rather than 'maxFactor'"
124  << exit(FatalIOError);
125  }
126 
127  const scalar pMaxFactor(dict.lookup<scalar>("maxFactor"));
128  pMax_.value() = pMaxFactor*pMax;
129  limitMaxP_ = true;
130  }
131  }
132 
133  if (limitMinP_ || limitMaxP_)
134  {
135  if (limitMinP_)
136  {
137  Info<< " min " << pMin_.value() << nl;
138  }
139 
140  if (limitMaxP_)
141  {
142  Info<< " max " << pMax_.value() << nl;
143  }
144 
145  Info << endl;
146  }
147 }
148 
149 
150 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
151 
153 (
154  const word& name,
155  const word& modelType,
156  const dictionary& dict,
157  const fvMesh& mesh
158 )
159 :
160  fvConstraint(name, modelType, dict, mesh),
161  pMin_("pMin", dimPressure, 0),
162  pMax_("pMax", dimPressure, great),
163  limitMinP_(false),
164  limitMaxP_(false)
165 {
166  readCoeffs();
167 }
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 {
174  return wordList(1, "p");
175 }
176 
177 
179 {
180  bool constrained = false;
181 
182  if (limitMinP_ || limitMaxP_)
183  {
184  if (limitMinP_)
185  {
186  const scalar pMin = min(p).value();
187 
188  if (pMin < pMin_.value())
189  {
190  Info<< "limitPressure: p min " << pMin << endl;
191  p = max(p, pMin_);
192  constrained = true;
193  }
194  }
195 
196  if (limitMaxP_)
197  {
198  const scalar pMax = max(p).value();
199 
200  if (pMax > pMax_.value())
201  {
202  Info<< "limitPressure: p max " << pMax << endl;
203  p = min(p, pMax_);
204  constrained = true;
205  }
206  }
207 
209  }
210 
211  return constrained;
212 }
213 
214 
216 {
217  if (fvConstraint::read(dict))
218  {
219  readCoeffs();
220  return true;
221  }
222  else
223  {
224  return false;
225  }
226 }
227 
228 
229 // ************************************************************************* //
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#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:156
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionSet dimPressure
dimensionedScalar pMin("pMin", dimPressure, mixture)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
virtual bool constrain(volScalarField &he) const
Constrain the energy field.
A class for handling words, derived from string.
Definition: word.H:59
labelList fv(nPoints)
static const char nl
Definition: Ostream.H:260
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
addToRunTimeSelectionTable(fvConstraint, fixedTemperatureConstraint, dictionary)
List< word > wordList
A List of words.
Definition: fileName.H:54
label patchi
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
limitPressure(const word &name, const word &modelType, const dictionary &dict, const fvMesh &mesh)
Construct from components.
void correctBoundaryConditions()
Correct boundary field.
messageStream Info
virtual bool read(const dictionary &dict)
Read dictionary.
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:143
Finite volume options abstract base class.
Definition: fvConstraint.H:54
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
IOerror FatalIOError