uniformTotalPressureFvPatchScalarField.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-2021 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// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36
37Foam::uniformTotalPressureFvPatchScalarField::
38uniformTotalPressureFvPatchScalarField
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_()
51{}
52
53
54Foam::uniformTotalPressureFvPatchScalarField::
55uniformTotalPressureFvPatchScalarField
56(
57 const fvPatch& p,
59 const dictionary& dict
60)
61:
62 fixedValueFvPatchScalarField(p, iF, dict, false),
63 UName_(dict.getOrDefault<word>("U", "U")),
64 phiName_(dict.getOrDefault<word>("phi", "phi")),
65 rhoName_(dict.getOrDefault<word>("rho", "rho")),
66 psiName_(dict.getOrDefault<word>("psi", "none")),
67 gamma_(psiName_ != "none" ? dict.get<scalar>("gamma") : 1),
68 p0_(Function1<scalar>::New("p0", dict, &db()))
69{
70 if (dict.found("value"))
71 {
73 (
74 scalarField("value", dict, p.size())
75 );
76 }
77 else
78 {
79 const scalar t = this->db().time().timeOutputValue();
80 fvPatchScalarField::operator==(p0_->value(t));
81 }
82}
83
84
85Foam::uniformTotalPressureFvPatchScalarField::
86uniformTotalPressureFvPatchScalarField
87(
89 const fvPatch& p,
91 const fvPatchFieldMapper& mapper
92)
93:
94 fixedValueFvPatchScalarField(p, iF), // Don't map
95 UName_(ptf.UName_),
96 phiName_(ptf.phiName_),
97 rhoName_(ptf.rhoName_),
98 psiName_(ptf.psiName_),
99 gamma_(ptf.gamma_),
100 p0_(ptf.p0_.clone())
101{
102 patchType() = ptf.patchType();
103
104 // Set the patch pressure to the current total pressure
105 // This is not ideal but avoids problems with the creation of patch faces
106 const scalar t = this->db().time().timeOutputValue();
107 fvPatchScalarField::operator==(p0_->value(t));
108}
109
110
111Foam::uniformTotalPressureFvPatchScalarField::
112uniformTotalPressureFvPatchScalarField
113(
115)
116:
117 fixedValueFvPatchScalarField(ptf),
118 UName_(ptf.UName_),
119 phiName_(ptf.phiName_),
120 rhoName_(ptf.rhoName_),
121 psiName_(ptf.psiName_),
122 gamma_(ptf.gamma_),
123 p0_(ptf.p0_.clone())
124{}
125
126
127Foam::uniformTotalPressureFvPatchScalarField::
128uniformTotalPressureFvPatchScalarField
129(
132)
133:
134 fixedValueFvPatchScalarField(ptf, iF),
135 UName_(ptf.UName_),
136 phiName_(ptf.phiName_),
137 rhoName_(ptf.rhoName_),
138 psiName_(ptf.psiName_),
139 gamma_(ptf.gamma_),
140 p0_(ptf.p0_.clone())
141{}
142
143
144// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
145
147(
148 const vectorField& Up
149)
150{
151 if (updated())
152 {
153 return;
154 }
155
156 scalar p0 = p0_->value(this->db().time().timeOutputValue());
157
158 const fvsPatchField<scalar>& phip =
159 patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
160
161 if (internalField().dimensions() == dimPressure)
162 {
163 if (psiName_ == "none")
164 {
165 // Variable density and low-speed compressible flow
166
168 patch().lookupPatchField<volScalarField, scalar>(rhoName_);
169
170 operator==(p0 - 0.5*rho*(1.0 - pos0(phip))*magSqr(Up));
171 }
172 else
173 {
174 // High-speed compressible flow
175
176 const fvPatchField<scalar>& psip =
177 patch().lookupPatchField<volScalarField, scalar>(psiName_);
178
179 if (gamma_ > 1)
180 {
181 scalar gM1ByG = (gamma_ - 1)/gamma_;
182
183 operator==
184 (
185 p0
186 /pow
187 (
188 (1.0 + 0.5*psip*gM1ByG*(1.0 - pos0(phip))*magSqr(Up)),
189 1.0/gM1ByG
190 )
191 );
192 }
193 else
194 {
195 operator==(p0/(1.0 + 0.5*psip*(1.0 - pos0(phip))*magSqr(Up)));
196 }
197 }
198
199 }
200 else if (internalField().dimensions() == dimPressure/dimDensity)
201 {
202 // Incompressible flow
203 operator==(p0 - 0.5*(1.0 - pos0(phip))*magSqr(Up));
204 }
205 else
206 {
208 << " Incorrect pressure dimensions " << internalField().dimensions()
209 << nl
210 << " Should be " << dimPressure
211 << " for compressible/variable density flow" << nl
212 << " or " << dimPressure/dimDensity
213 << " for incompressible flow," << nl
214 << " on patch " << this->patch().name()
215 << " of field " << this->internalField().name()
216 << " in file " << this->internalField().objectPath()
217 << exit(FatalError);
218 }
219
220 fixedValueFvPatchScalarField::updateCoeffs();
221}
222
223
225{
226 updateCoeffs(patch().lookupPatchField<volVectorField, vector>(UName_));
227}
228
229
231{
233 os.writeEntryIfDifferent<word>("U", "U", UName_);
234 os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
235 os.writeEntry("rho", rhoName_);
236 os.writeEntry("psi", psiName_);
237 os.writeEntry("gamma", gamma_);
238 p0_->writeData(os);
239 writeEntry("value", os);
240}
241
242
243// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
244
245namespace Foam
246{
248 (
251 );
252}
253
254// ************************************************************************* //
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...
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:96
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
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
friend bool operator==(const refineCell &rc1, const refineCell &rc2)
Definition: refineCell.H:97
This boundary condition provides a time-varying form of the uniform total pressure boundary condition...
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 > &)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
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.