totalPressureFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
31 #include "fvPatchFieldMapper.H"
32 #include "volFields.H"
33 #include "surfaceFields.H"
34 
35 
36 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
37 
39 (
40  const fvPatch& p,
42 )
43 :
44  fixedValueFvPatchScalarField(p, iF),
45  UName_("U"),
46  phiName_("phi"),
47  rhoName_("rho"),
48  psiName_("none"),
49  gamma_(0.0),
50  p0_(p.size(), Zero)
51 {}
52 
53 
55 (
56  const fvPatch& p,
58  const dictionary& dict
59 )
60 :
61  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
62  UName_(dict.getOrDefault<word>("U", "U")),
63  phiName_(dict.getOrDefault<word>("phi", "phi")),
64  rhoName_(dict.getOrDefault<word>("rho", "rho")),
65  psiName_(dict.getOrDefault<word>("psi", "none")),
66  gamma_(psiName_ != "none" ? dict.get<scalar>("gamma") : 1),
67  p0_("p0", dict, p.size())
68 {
69  if (!this->readValueEntry(dict))
70  {
72  }
73 }
74 
75 
77 (
78  const totalPressureFvPatchScalarField& ptf,
79  const fvPatch& p,
80  const DimensionedField<scalar, volMesh>& iF,
81  const fvPatchFieldMapper& mapper
82 )
83 :
84  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
85  UName_(ptf.UName_),
86  phiName_(ptf.phiName_),
87  rhoName_(ptf.rhoName_),
88  psiName_(ptf.psiName_),
89  gamma_(ptf.gamma_),
90  p0_(ptf.p0_, mapper)
91 {}
92 
93 
95 (
97 )
98 :
99  fixedValueFvPatchScalarField(tppsf),
100  UName_(tppsf.UName_),
101  phiName_(tppsf.phiName_),
102  rhoName_(tppsf.rhoName_),
103  psiName_(tppsf.psiName_),
104  gamma_(tppsf.gamma_),
105  p0_(tppsf.p0_)
106 {}
107 
108 
110 (
111  const totalPressureFvPatchScalarField& tppsf,
113 )
114 :
115  fixedValueFvPatchScalarField(tppsf, iF),
116  UName_(tppsf.UName_),
117  phiName_(tppsf.phiName_),
118  rhoName_(tppsf.rhoName_),
119  psiName_(tppsf.psiName_),
120  gamma_(tppsf.gamma_),
121  p0_(tppsf.p0_)
122 {}
123 
124 
125 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
126 
128 (
129  const fvPatchFieldMapper& m
130 )
131 {
132  fixedValueFvPatchScalarField::autoMap(m);
133  p0_.autoMap(m);
134 }
135 
136 
138 (
139  const fvPatchScalarField& ptf,
140  const labelList& addr
141 )
142 {
143  fixedValueFvPatchScalarField::rmap(ptf, addr);
144 
145  const totalPressureFvPatchScalarField& tiptf =
146  refCast<const totalPressureFvPatchScalarField>(ptf);
147 
148  p0_.rmap(tiptf.p0_, addr);
149 }
150 
151 
153 (
154  const scalarField& p0p,
155  const vectorField& Up
156 )
157 {
158  if (updated())
159  {
160  return;
161  }
162 
163  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
164 
165  if (internalField().dimensions() == dimPressure)
166  {
167  if (psiName_ == "none")
168  {
169  // Variable density and low-speed compressible flow
170 
171  const auto& rho =
172  patch().lookupPatchField<volScalarField>(rhoName_);
173 
174  operator==(p0p - 0.5*rho*(neg(phip))*magSqr(Up));
175  }
176  else
177  {
178  // High-speed compressible flow
179 
180  const auto& psip =
181  patch().lookupPatchField<volScalarField>(psiName_);
182 
183  if (gamma_ > 1)
184  {
185  scalar gM1ByG = (gamma_ - 1)/gamma_;
186 
187  operator==
188  (
189  p0p
190  /pow
191  (
192  (1.0 + 0.5*psip*gM1ByG*(neg(phip))*magSqr(Up)),
193  1.0/gM1ByG
194  )
195  );
196  }
197  else
198  {
199  operator==(p0p/(1.0 + 0.5*psip*(neg(phip))*magSqr(Up)));
200  }
201  }
202 
203  }
204  else if (internalField().dimensions() == dimPressure/dimDensity)
205  {
206  // Incompressible flow
207  operator==(p0p - 0.5*(neg(phip))*magSqr(Up));
208  }
209  else
210  {
212  << " Incorrect pressure dimensions " << internalField().dimensions()
213  << nl
214  << " Should be " << dimPressure
215  << " for compressible/variable density flow" << nl
216  << " or " << dimPressure/dimDensity
217  << " for incompressible flow," << nl
218  << " on patch " << this->patch().name()
219  << " of field " << this->internalField().name()
220  << " in file " << this->internalField().objectPath()
222  }
223 
224  fixedValueFvPatchScalarField::updateCoeffs();
225 }
226 
227 
229 {
230  updateCoeffs
231  (
232  p0(),
233  patch().lookupPatchField<volVectorField>(UName())
234  );
235 }
236 
237 
239 {
241  os.writeEntryIfDifferent<word>("U", "U", UName_);
242  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
243  os.writeEntry("rho", rhoName_);
244  os.writeEntry("psi", psiName_);
245  os.writeEntry("gamma", gamma_);
246  p0_.writeEntry("p0", os);
248 }
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 namespace Foam
254 {
256  (
258  totalPressureFvPatchScalarField
259  );
260 }
261 
262 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
const word UName(propsDict.getOrDefault< word >("U", "U"))
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
fvPatchField< scalar > fvPatchScalarField
A class for handling words, derived from Foam::string.
Definition: word.H:63
A FieldMapper for finite-volume patch fields.
const dimensionSet dimPressure
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
OBJstream os(runTime.globalPath()/outputName)
const dimensionSet dimDensity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
totalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:391
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
This boundary condition provides a total pressure condition. Four variants are possible: ...