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-2024 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 {
38  (
42  );
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::fv::limitPressure::readCoeffs(const dictionary& dict)
50 {
51  pName_ = dict.lookupOrDefault<word>("p", "p");
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>(pName_);
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
73  (
74  pbf[patchi].fixesValue()
75  || isA<calculatedFvPatchField<scalar>>(pbf[patchi])
76  )
77  {
78  pLimits = true;
79 
80  pMin = min(pMin, min(pbf[patchi]));
81  pMax = max(pMax, max(pbf[patchi]));
82  }
83  }
84 
85  reduce(pLimits, andOp<bool>());
86  if (pLimits)
87  {
88  reduce(pMax, maxOp<scalar>());
89  reduce(pMin, minOp<scalar>());
90  }
91 
92  if (dict.found("min"))
93  {
94  pMin_.value() = dict.lookup<scalar>("min");
95  limitMinP_ = true;
96  }
97  else if (dict.found("minFactor"))
98  {
99  if (!pLimits)
100  {
102  << "'minFactor' specified rather than 'min'" << nl
103  << " but the corresponding reference pressure cannot"
104  " be evaluated from the boundary conditions." << nl
105  << " Please specify 'min' rather than 'minFactor'"
106  << exit(FatalIOError);
107  }
108 
109  const scalar pMinFactor(dict.lookup<scalar>("minFactor"));
110  pMin_.value() = pMinFactor*pMin;
111  limitMinP_ = true;
112  }
113 
114  if (dict.found("max"))
115  {
116  pMax_.value() = dict.lookup<scalar>("max");
117  limitMaxP_ = true;
118  }
119  else if (dict.found("maxFactor"))
120  {
121  if (!pLimits)
122  {
124  << "'maxFactor' specified rather than 'max'" << nl
125  << " but the corresponding reference pressure cannot"
126  " be evaluated from the boundary conditions." << nl
127  << " Please specify 'max' rather than 'maxFactor'"
128  << exit(FatalIOError);
129  }
130 
131  const scalar pMaxFactor(dict.lookup<scalar>("maxFactor"));
132  pMax_.value() = pMaxFactor*pMax;
133  limitMaxP_ = true;
134  }
135  }
136 
137  if (limitMinP_ || limitMaxP_)
138  {
139  if (limitMinP_)
140  {
141  Info<< " min " << pMin_.value() << nl;
142  }
143 
144  if (limitMaxP_)
145  {
146  Info<< " max " << pMax_.value() << nl;
147  }
148 
149  Info << endl;
150  }
151 }
152 
153 
154 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
155 
157 (
158  const word& name,
159  const word& modelType,
160  const fvMesh& mesh,
161  const dictionary& dict
162 )
163 :
164  fvConstraint(name, modelType, mesh, dict),
165  pName_(word::null),
166  pMin_("pMin", dimPressure, 0),
167  pMax_("pMax", dimPressure, great),
168  limitMinP_(false),
169  limitMaxP_(false)
170 {
171  readCoeffs(coeffs(dict));
172 }
173 
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 {
179  return wordList(1, pName_);
180 }
181 
182 
184 {
185  bool constrained = false;
186 
187  if (limitMinP_ || limitMaxP_)
188  {
189  if (limitMinP_)
190  {
191  const scalar pMin = min(p).value();
192 
193  if (pMin < pMin_.value())
194  {
195  Info<< "limitPressure: min " << pMin << endl;
196  p = max(p, pMin_);
197  constrained = true;
198  }
199  }
200 
201  if (limitMaxP_)
202  {
203  const scalar pMax = max(p).value();
204 
205  if (pMax > pMax_.value())
206  {
207  Info<< "limitPressure: max " << pMax << endl;
208  p = min(p, pMax_);
209  constrained = true;
210  }
211  }
212 
213  p.correctBoundaryConditions();
214  }
215 
216  return constrained;
217 }
218 
219 
221 {
222  return true;
223 }
224 
225 
227 {}
228 
229 
231 {}
232 
233 
235 {}
236 
237 
239 {
241  {
242  readCoeffs(coeffs(dict));
243  return true;
244  }
245  else
246  {
247  return false;
248  }
249 }
250 
251 
252 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
GeometricBoundaryField< Type, GeoMesh, PrimitiveField > Boundary
Type of the boundary field.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const Type & value() const
Return const reference to value.
Finite volume options abstract base class.
Definition: fvConstraint.H:57
virtual bool read(const dictionary &dict)
Read source dictionary.
Definition: fvConstraint.C:164
const dictionary & coeffs(const dictionary &) const
Return the coefficients sub-dictionary.
Definition: fvConstraintI.H:41
const fvMesh & mesh() const
Return const access to the mesh database.
Definition: fvConstraintI.H:34
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Limits the specified pressure field to be between specified minimum and maximum limits.
Definition: limitPressure.H:69
virtual bool movePoints()
Update for mesh motion.
virtual bool constrain(volScalarField &p) const
Constrain the pressure field.
limitPressure(const word &name, const word &modelType, const fvMesh &mesh, const dictionary &dict)
Construct from components.
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
virtual void distribute(const polyDistributionMap &)
Redistribute or update using the given distribution map.
virtual bool read(const dictionary &dict)
Read dictionary.
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
virtual wordList constrainedFields() const
Return the list of fields constrained by the fvConstraint.
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type and name.
Class containing mesh-to-mesh mapping information after a mesh distribution where we send parts of me...
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:51
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
label patchi
addToRunTimeSelectionTable(fvConstraint, bound, dictionary)
defineTypeNameAndDebug(bound, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< word > wordList
A List of words.
Definition: fileName.H:54
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:171
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
static const char nl
Definition: Ostream.H:267
labelList fv(nPoints)
dictionary dict
volScalarField & p