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 -------------------------------------------------------------------------------
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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
36 
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 
56 (
57  const fvPatch& p,
59  const dictionary& dict
60 )
61 :
62  fixedValueFvPatchScalarField(p, iF, dict, IOobjectOption::NO_READ),
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 (!this->readValueEntry(dict))
71  {
72  const scalar t = this->db().time().timeOutputValue();
74  }
75 }
76 
77 
80 (
81  const uniformTotalPressureFvPatchScalarField& ptf,
82  const fvPatch& p,
83  const DimensionedField<scalar, volMesh>& iF,
84  const fvPatchFieldMapper& mapper
85 )
86 :
87  fixedValueFvPatchScalarField(p, iF), // Don't map
88  UName_(ptf.UName_),
89  phiName_(ptf.phiName_),
90  rhoName_(ptf.rhoName_),
91  psiName_(ptf.psiName_),
92  gamma_(ptf.gamma_),
93  p0_(ptf.p0_.clone())
94 {
95  patchType() = ptf.patchType();
96 
97  // Set the patch pressure to the current total pressure
98  // This is not ideal but avoids problems with the creation of patch faces
99  const scalar t = this->db().time().timeOutputValue();
100  fvPatchScalarField::operator==(p0_->value(t));
101 }
102 
103 
106 (
108 )
109 :
110  fixedValueFvPatchScalarField(ptf),
111  UName_(ptf.UName_),
112  phiName_(ptf.phiName_),
113  rhoName_(ptf.rhoName_),
114  psiName_(ptf.psiName_),
115  gamma_(ptf.gamma_),
116  p0_(ptf.p0_.clone())
117 {}
118 
119 
122 (
125 )
126 :
127  fixedValueFvPatchScalarField(ptf, iF),
128  UName_(ptf.UName_),
129  phiName_(ptf.phiName_),
130  rhoName_(ptf.rhoName_),
131  psiName_(ptf.psiName_),
132  gamma_(ptf.gamma_),
133  p0_(ptf.p0_.clone())
134 {}
135 
136 
137 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
138 
140 (
141  const vectorField& Up
142 )
143 {
144  if (updated())
145  {
146  return;
147  }
148 
149  scalar p0 = p0_->value(this->db().time().timeOutputValue());
150 
151  const auto& phip = patch().lookupPatchField<surfaceScalarField>(phiName_);
152 
153  if (internalField().dimensions() == dimPressure)
154  {
155  if (psiName_ == "none")
156  {
157  // Variable density and low-speed compressible flow
158 
159  const auto& rho =
160  patch().lookupPatchField<volScalarField>(rhoName_);
161 
162  operator==(p0 - 0.5*rho*(neg(phip))*magSqr(Up));
163  }
164  else
165  {
166  // High-speed compressible flow
167 
168  const auto& psip =
169  patch().lookupPatchField<volScalarField>(psiName_);
170 
171  if (gamma_ > 1)
172  {
173  scalar gM1ByG = (gamma_ - 1)/gamma_;
174 
175  operator==
176  (
177  p0
178  /pow
179  (
180  (1.0 + 0.5*psip*gM1ByG*(neg(phip))*magSqr(Up)),
181  1.0/gM1ByG
182  )
183  );
184  }
185  else
186  {
187  operator==(p0/(1.0 + 0.5*psip*(neg(phip))*magSqr(Up)));
188  }
189  }
190 
191  }
192  else if (internalField().dimensions() == dimPressure/dimDensity)
193  {
194  // Incompressible flow
195  operator==(p0 - 0.5*(neg(phip))*magSqr(Up));
196  }
197  else
198  {
200  << " Incorrect pressure dimensions " << internalField().dimensions()
201  << nl
202  << " Should be " << dimPressure
203  << " for compressible/variable density flow" << nl
204  << " or " << dimPressure/dimDensity
205  << " for incompressible flow," << nl
206  << " on patch " << this->patch().name()
207  << " of field " << this->internalField().name()
208  << " in file " << this->internalField().objectPath()
210  }
211 
212  fixedValueFvPatchScalarField::updateCoeffs();
213 }
214 
217 {
218  updateCoeffs(patch().lookupPatchField<volVectorField>(UName_));
219 }
220 
221 
223 {
225  os.writeEntryIfDifferent<word>("U", "U", UName_);
226  os.writeEntryIfDifferent<word>("phi", "phi", phiName_);
227  os.writeEntry("rho", rhoName_);
228  os.writeEntry("psi", psiName_);
229  os.writeEntry("gamma", gamma_);
230  p0_->writeData(os);
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
236 
237 namespace Foam
238 {
240  (
242  uniformTotalPressureFvPatchScalarField
243  );
244 }
245 
246 // ************************************************************************* //
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.
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
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
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
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
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.
virtual void operator==(const fvPatchField< scalar > &)
Definition: fvPatchField.C:546
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
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
uniformTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
OBJstream os(runTime.globalPath()/outputName)
This boundary condition provides a time-varying form of the uniform total pressure boundary condition...
const dimensionSet dimDensity
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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.