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-2022 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  pName_ = dict.lookupOrDefault<word>("p", "p");
54 
55  if (dict.found("min") && dict.found("max"))
56  {
57  pMin_.value() = dict.lookup<scalar>("min");
58  limitMinP_ = true;
59 
60  pMax_.value() = dict.lookup<scalar>("max");
61  limitMaxP_ = true;
62  }
63  else
64  {
65  const volScalarField& p = mesh().lookupObject<volScalarField>(pName_);
66  const volScalarField::Boundary& pbf = p.boundaryField();
67 
68  bool pLimits = false;
69  scalar pMin = vGreat;
70  scalar pMax = -vGreat;
71 
72  forAll(pbf, patchi)
73  {
74  if
75  (
76  pbf[patchi].fixesValue()
77  || isA<calculatedFvPatchField<scalar>>(pbf[patchi])
78  )
79  {
80  pLimits = true;
81 
82  pMin = min(pMin, min(pbf[patchi]));
83  pMax = max(pMax, max(pbf[patchi]));
84  }
85  }
86 
87  reduce(pLimits, andOp<bool>());
88  if (pLimits)
89  {
90  reduce(pMax, maxOp<scalar>());
91  reduce(pMin, minOp<scalar>());
92  }
93 
94  if (dict.found("min"))
95  {
96  pMin_.value() = dict.lookup<scalar>("min");
97  limitMinP_ = true;
98  }
99  else if (dict.found("minFactor"))
100  {
101  if (!pLimits)
102  {
104  << "'minFactor' specified rather than 'min'" << nl
105  << " but the corresponding reference pressure cannot"
106  " be evaluated from the boundary conditions." << nl
107  << " Please specify 'min' rather than 'minFactor'"
108  << exit(FatalIOError);
109  }
110 
111  const scalar pMinFactor(dict.lookup<scalar>("minFactor"));
112  pMin_.value() = pMinFactor*pMin;
113  limitMinP_ = true;
114  }
115 
116  if (dict.found("max"))
117  {
118  pMax_.value() = dict.lookup<scalar>("max");
119  limitMaxP_ = true;
120  }
121  else if (dict.found("maxFactor"))
122  {
123  if (!pLimits)
124  {
126  << "'maxFactor' specified rather than 'max'" << nl
127  << " but the corresponding reference pressure cannot"
128  " be evaluated from the boundary conditions." << nl
129  << " Please specify 'max' rather than 'maxFactor'"
130  << exit(FatalIOError);
131  }
132 
133  const scalar pMaxFactor(dict.lookup<scalar>("maxFactor"));
134  pMax_.value() = pMaxFactor*pMax;
135  limitMaxP_ = true;
136  }
137  }
138 
139  if (limitMinP_ || limitMaxP_)
140  {
141  if (limitMinP_)
142  {
143  Info<< " min " << pMin_.value() << nl;
144  }
145 
146  if (limitMaxP_)
147  {
148  Info<< " max " << pMax_.value() << nl;
149  }
150 
151  Info << endl;
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
157 
159 (
160  const word& name,
161  const word& modelType,
162  const dictionary& dict,
163  const fvMesh& mesh
164 )
165 :
166  fvConstraint(name, modelType, dict, mesh),
167  pName_(word::null),
168  pMin_("pMin", dimPressure, 0),
169  pMax_("pMax", dimPressure, great),
170  limitMinP_(false),
171  limitMaxP_(false)
172 {
173  readCoeffs();
174 }
175 
176 
177 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
178 
180 {
181  return wordList(1, pName_);
182 }
183 
184 
186 {
187  bool constrained = false;
188 
189  if (limitMinP_ || limitMaxP_)
190  {
191  if (limitMinP_)
192  {
193  const scalar pMin = min(p).value();
194 
195  if (pMin < pMin_.value())
196  {
197  Info<< "limitPressure: min " << pMin << endl;
198  p = max(p, pMin_);
199  constrained = true;
200  }
201  }
202 
203  if (limitMaxP_)
204  {
205  const scalar pMax = max(p).value();
206 
207  if (pMax > pMax_.value())
208  {
209  Info<< "limitPressure: max " << pMax << endl;
210  p = min(p, pMax_);
211  constrained = true;
212  }
213  }
214 
216  }
217 
218  return constrained;
219 }
220 
221 
223 {
224  return true;
225 }
226 
227 
229 {}
230 
231 
233 {}
234 
235 
237 {}
238 
239 
240 
242 {
243  if (fvConstraint::read(dict))
244  {
245  readCoeffs();
246  return true;
247  }
248  else
249  {
250  return false;
251  }
252 }
253 
254 
255 // ************************************************************************* //
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:134
dictionary dict
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
defineTypeNameAndDebug(fixedTemperatureConstraint, 0)
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
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
const dimensionSet dimPressure
dimensionedScalar pMin("pMin", dimPressure, mixture)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of the boundary field.
fvMesh & mesh
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
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)
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
static const word null
An empty word.
Definition: word.H:77
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const bool constrained
Definition: pEqn.H:108
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
static const char nl
Definition: Ostream.H:260
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
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
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
virtual bool movePoints()
Update for mesh motion.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
Finite volume options abstract base class.
Definition: fvConstraint.H:56
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864
IOerror FatalIOError