pressure.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) 2012-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 "pressure.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace functionObjects
35 {
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
43 
44 Foam::word Foam::functionObjects::pressure::resultName() const
45 {
46  word rName;
47 
48  if (calcTotal_)
49  {
50  rName = "total(" + fieldName_ + ")";
51  }
52  else
53  {
54  rName = "static(" + fieldName_ + ")";
55  }
56 
57  if (calcCoeff_)
58  {
59  rName += "_coeff";
60  }
61 
62  return rName;
63 }
64 
65 
66 Foam::tmp<Foam::volScalarField> Foam::functionObjects::pressure::rhoScale
67 (
68  const volScalarField& p
69 ) const
70 {
71  if (p.dimensions() == dimPressure)
72  {
73  return p;
74  }
75  else
76  {
77  return dimensionedScalar("rhoInf", dimDensity, rhoInf_)*p;
78  }
79 }
80 
81 
82 Foam::tmp<Foam::volScalarField> Foam::functionObjects::pressure::rhoScale
83 (
84  const volScalarField& p,
85  const tmp<volScalarField>& tsf
86 ) const
87 {
88  if (p.dimensions() == dimPressure)
89  {
90  return lookupObject<volScalarField>(rhoName_)*tsf;
91  }
92  else
93  {
94  return dimensionedScalar("rhoInf", dimDensity, rhoInf_)*tsf;
95  }
96 }
97 
98 
99 Foam::tmp<Foam::volScalarField> Foam::functionObjects::pressure::pRef
100 (
101  const tmp<volScalarField>& tp
102 ) const
103 {
104  if (calcTotal_)
105  {
106  return tp + dimensionedScalar("pRef", dimPressure, pRef_);
107  }
108  else
109  {
110  return move(tp);
111  }
112 }
113 
114 
115 Foam::tmp<Foam::volScalarField> Foam::functionObjects::pressure::pDyn
116 (
117  const volScalarField& p,
118  const tmp<volScalarField>& tp
119 ) const
120 {
121  if (calcTotal_)
122  {
123  return
124  tp + rhoScale(p, 0.5*magSqr(lookupObject<volVectorField>(UName_)));
125  }
126  else
127  {
128  return move(tp);
129  }
130 }
131 
132 
134 Foam::functionObjects::pressure::coeff
135 (
136  const tmp<volScalarField>& tp
137 ) const
138 {
139  if (calcCoeff_)
140  {
141  tmp<volScalarField> tpCoeff(tp.ptr());
142  volScalarField& pCoeff = tpCoeff.ref();
143 
144  pCoeff -= dimensionedScalar("pInf", dimPressure, pInf_);
145 
146  const dimensionedScalar pSmall("pSmall", dimPressure, small);
147  const dimensionedVector U("U", dimVelocity, UInf_);
148  const dimensionedScalar rho("rho", dimDensity, rhoInf_);
149 
150  pCoeff /= 0.5*rho*magSqr(U) + pSmall;
151 
152  return tpCoeff;
153  }
154  else
155  {
156  return move(tp);
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
162 
163 bool Foam::functionObjects::pressure::calc()
164 {
165  if (foundObject<volScalarField>(fieldName_))
166  {
167  const volScalarField& p = lookupObject<volScalarField>(fieldName_);
168 
169  store
170  (
171  resultName_,
172  volScalarField::New(resultName_, coeff(pRef(pDyn(p, rhoScale(p)))))
173  );
174 
175  return true;
176  }
177  else
178  {
179  cannotFindObject<volScalarField>(fieldName_);
180 
181  return false;
182  }
183 }
184 
185 
186 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
187 
189 (
190  const word& name,
191  const Time& runTime,
192  const dictionary& dict
193 )
194 :
195  fieldExpression(name, runTime, dict, typeName, "p"),
196  UName_("U"),
197  rhoName_("rho"),
198  calcTotal_(false),
199  pRef_(0),
200  calcCoeff_(false),
201  pInf_(0),
202  UInf_(Zero),
203  rhoInf_(1)
204 {
205  read(dict);
206 
207  dimensionSet pDims(dimPressure);
208 
209  if (calcCoeff_)
210  {
211  pDims /= dimPressure;
212  }
213 }
214 
215 
216 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
217 
219 {}
220 
221 
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223 
225 {
226  dict.readIfPresent("U", UName_);
227  dict.readIfPresent("rho", rhoName_);
228 
229  if (rhoName_ == "rhoInf")
230  {
231  dict.lookup("rhoInf") >> rhoInf_;
232  }
233 
234  dict.lookup("calcTotal") >> calcTotal_;
235  if (calcTotal_)
236  {
237  dict.lookup("pRef") >> pRef_;
238  }
239 
240  dict.lookup("calcCoeff") >> calcCoeff_;
241  if (calcCoeff_)
242  {
243  dict.lookup("pInf") >> pInf_;
244  dict.lookup("UInf") >> UInf_;
245  dict.lookup("rhoInf") >> rhoInf_;
246 
247  scalar zeroCheck = 0.5*rhoInf_*magSqr(UInf_) + pInf_;
248 
249  if (mag(zeroCheck) < rootVSmall)
250  {
252  << type() << " " << name() << ": "
253  << "Coefficient calculation requested, but reference "
254  << "pressure level is zero. Please check the supplied "
255  << "values of pInf, UInf and rhoInf" << endl;
256  }
257  }
258 
259  resultName_ = dict.lookupOrDefault<word>("result", resultName());
260 
261  return true;
262 }
263 
264 
265 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
Dimension set for the base types.
Definition: dimensionSet.H:125
Abstract base-class for Time/database functionObjects.
const word fieldName_
Name of field to process.
Includes tools to manipulate the pressure into different forms.
Definition: pressure.H:257
pressure(const word &name, const Time &runTime, const dictionary &)
Construct from Time and dictionary.
Definition: pressure.C:189
virtual ~pressure()
Destructor.
Definition: pressure.C:218
virtual bool read(const dictionary &)
Read the pressure data.
Definition: pressure.C:224
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
addToRunTimeSelectionTable(functionObject, adjustTimeStepToCombustion, dictionary)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimVelocity
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
dictionary dict
volScalarField & p