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-------------------------------------------------------------------------------
11License
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, false),
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 (dict.found("value"))
70 {
72 (
73 scalarField("value", dict, p.size())
74 );
75 }
76 else
77 {
79 }
80}
81
82
84(
86 const fvPatch& p,
88 const fvPatchFieldMapper& mapper
89)
90:
91 fixedValueFvPatchScalarField(ptf, p, iF, mapper),
92 UName_(ptf.UName_),
93 phiName_(ptf.phiName_),
94 rhoName_(ptf.rhoName_),
95 psiName_(ptf.psiName_),
96 gamma_(ptf.gamma_),
97 p0_(ptf.p0_, mapper)
98{}
99
100
102(
104)
105:
106 fixedValueFvPatchScalarField(tppsf),
107 UName_(tppsf.UName_),
108 phiName_(tppsf.phiName_),
109 rhoName_(tppsf.rhoName_),
110 psiName_(tppsf.psiName_),
111 gamma_(tppsf.gamma_),
112 p0_(tppsf.p0_)
113{}
114
115
117(
120)
121:
122 fixedValueFvPatchScalarField(tppsf, iF),
123 UName_(tppsf.UName_),
124 phiName_(tppsf.phiName_),
125 rhoName_(tppsf.rhoName_),
126 psiName_(tppsf.psiName_),
127 gamma_(tppsf.gamma_),
128 p0_(tppsf.p0_)
129{}
130
131
132// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
133
135(
136 const fvPatchFieldMapper& m
137)
138{
139 fixedValueFvPatchScalarField::autoMap(m);
140 p0_.autoMap(m);
141}
142
143
145(
146 const fvPatchScalarField& ptf,
147 const labelList& addr
148)
149{
150 fixedValueFvPatchScalarField::rmap(ptf, addr);
151
153 refCast<const totalPressureFvPatchScalarField>(ptf);
154
155 p0_.rmap(tiptf.p0_, addr);
156}
157
158
160(
161 const scalarField& p0p,
162 const vectorField& Up
163)
164{
165 if (updated())
166 {
167 return;
168 }
169
170 const fvsPatchField<scalar>& phip =
171 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
172
173 if (internalField().dimensions() == dimPressure)
174 {
175 if (psiName_ == "none")
176 {
177 // Variable density and low-speed compressible flow
178
180 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
181
182 operator==(p0p - 0.5*rho*(1.0 - pos0(phip))*magSqr(Up));
183 }
184 else
185 {
186 // High-speed compressible flow
187
188 const fvPatchField<scalar>& psip =
189 patch().lookupPatchField<volScalarField, scalar>(psiName_);
190
191 if (gamma_ > 1)
192 {
193 scalar gM1ByG = (gamma_ - 1)/gamma_;
194
195 operator==
196 (
197 p0p
198 /pow
199 (
200 (1.0 + 0.5*psip*gM1ByG*(1.0 - pos0(phip))*magSqr(Up)),
201 1.0/gM1ByG
202 )
203 );
204 }
205 else
206 {
207 operator==(p0p/(1.0 + 0.5*psip*(1.0 - pos0(phip))*magSqr(Up)));
208 }
209 }
210
211 }
212 else if (internalField().dimensions() == dimPressure/dimDensity)
213 {
214 // Incompressible flow
215 operator==(p0p - 0.5*(1.0 - pos0(phip))*magSqr(Up));
216 }
217 else
218 {
220 << " Incorrect pressure dimensions " << internalField().dimensions()
221 << nl
222 << " Should be " << dimPressure
223 << " for compressible/variable density flow" << nl
224 << " or " << dimPressure/dimDensity
225 << " for incompressible flow," << nl
226 << " on patch " << this->patch().name()
227 << " of field " << this->internalField().name()
228 << " in file " << this->internalField().objectPath()
229 << exit(FatalError);
230 }
231
232 fixedValueFvPatchScalarField::updateCoeffs();
233}
234
235
237{
238 updateCoeffs
239 (
240 p0(),
241 patch().lookupPatchField<volVectorField, vector>(UName())
242 );
243}
244
245
247{
249 os.writeEntryIfDifferent<word>("U", "U", UName_);
250 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
251 os.writeEntry("rho", rhoName_);
252 os.writeEntry("psi", psiName_);
253 os.writeEntry("gamma", gamma_);
254 p0_.writeEntry("p0", os);
255 writeEntry("value", os);
256}
257
258
259// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260
261namespace Foam
262{
264 (
267 );
268}
269
270// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:239
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:251
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:408
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:71
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:79
This boundary condition provides a total pressure condition. Four variants are possible:
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition: word.H:68
volScalarField & p
const volScalarField & p0
Definition: EEqn.H:36
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
const dimensionSet dimPressure
dimensionedScalar pos0(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
error FatalError
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
dictionary dict
Foam::surfaceFields.
const word UName(propsDict.getOrDefault< word >("U", "U"))