uniformTotalPressureFvPatchScalarField.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) 2011-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 
28 #include "fieldMapper.H"
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  fixedValueFvPatchScalarField(p, iF, dict, false),
43  UName_(dict.lookupOrDefault<word>("U", "U")),
44  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
45  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
46  psiName_(dict.lookupOrDefault<word>("psi", "none")),
47  gamma_(psiName_ != "none" ? dict.lookup<scalar>("gamma", dimless) : 1),
48  p0_
49  (
50  Function1<scalar>::New
51  (
52  "p0",
53  db().time().userUnits(),
54  iF.dimensions(),
55  dict
56  )
57  )
58 {
59  if (dict.found("value"))
60  {
62  (
63  scalarField("value", iF.dimensions(), dict, p.size())
64  );
65  }
66  else
67  {
68  fvPatchScalarField::operator==(p0_->value(db().time().value()));
69  }
70 }
71 
72 
75 (
77  const fvPatch& p,
79  const fieldMapper& mapper
80 )
81 :
82  fixedValueFvPatchScalarField(ptf, p, iF, mapper, false), // Don't map
83  UName_(ptf.UName_),
84  phiName_(ptf.phiName_),
85  rhoName_(ptf.rhoName_),
86  psiName_(ptf.psiName_),
87  gamma_(ptf.gamma_),
88  p0_(ptf.p0_, false)
89 {
90  // Set the patch pressure to the current total pressure
91  // This is not ideal but avoids problems with the creation of patch faces
92  fvPatchScalarField::operator==(p0_->value(this->db().time().value()));
93 }
94 
95 
98 (
101 )
102 :
103  fixedValueFvPatchScalarField(ptf, iF),
104  UName_(ptf.UName_),
105  phiName_(ptf.phiName_),
106  rhoName_(ptf.rhoName_),
107  psiName_(ptf.psiName_),
108  gamma_(ptf.gamma_),
109  p0_(ptf.p0_, false)
110 {}
111 
112 
113 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 
116 (
117  const vectorField& Up
118 )
119 {
120  if (updated())
121  {
122  return;
123  }
124 
125  scalar p0 = p0_->value(this->db().time().value());
126 
127  const fvsPatchField<scalar>& phip =
128  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
129 
130  if (internalField().dimensions() == dimPressure)
131  {
132  if (psiName_ == "none")
133  {
134  // Variable density and low-speed compressible flow
135 
136  const fvPatchField<scalar>& rho =
137  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
138 
139  operator==(p0 - 0.5*rho*neg(phip)*magSqr(Up));
140  }
141  else
142  {
143  // High-speed compressible flow
144 
145  const fvPatchField<scalar>& psip =
146  patch().lookupPatchField<volScalarField, scalar>(psiName_);
147 
148  if (gamma_ > 1)
149  {
150  scalar gM1ByG = (gamma_ - 1)/gamma_;
151 
152  operator==
153  (
154  p0
155  /pow
156  (
157  (1.0 + 0.5*psip*gM1ByG*neg(phip)*magSqr(Up)),
158  1.0/gM1ByG
159  )
160  );
161  }
162  else
163  {
164  operator==(p0/(1.0 + 0.5*psip*neg(phip)*magSqr(Up)));
165  }
166  }
167 
168  }
169  else if (internalField().dimensions() == dimPressure/dimDensity)
170  {
171  // Incompressible flow
172  operator==(p0 - 0.5*neg(phip)*magSqr(Up));
173  }
174  else
175  {
177  << " Incorrect pressure dimensions " << internalField().dimensions()
178  << nl
179  << " Should be " << dimPressure
180  << " for compressible/variable density flow" << nl
181  << " or " << dimPressure/dimDensity
182  << " for incompressible flow," << nl
183  << " on patch " << this->patch().name()
184  << " of field " << this->internalField().name()
185  << " in file " << this->internalField().objectPath()
186  << exit(FatalError);
187  }
188 
189  fixedValueFvPatchScalarField::updateCoeffs();
190 }
191 
192 
194 {
195  updateCoeffs(patch().lookupPatchField<volVectorField, vector>(UName_));
196 }
197 
198 
200 {
202  writeEntryIfDifferent<word>(os, "U", "U", UName_);
203  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
204  writeEntry(os, "rho", rhoName_);
205  writeEntry(os, "psi", psiName_);
206  writeEntry(os, "gamma", gamma_);
207  writeEntry
208  (
209  os,
210  db().time().userUnits(),
211  internalField().dimensions(),
212  p0_()
213  );
214  writeEntry(os, "value", *this);
215 }
216 
217 
218 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
219 
220 namespace Foam
221 {
223  (
226  );
227 }
228 
229 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
Run-time selectable general function of one variable.
Definition: Function1.H:125
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:415
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
This boundary condition provides a time-varying form of the uniform total pressure boundary condition...
uniformTotalPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const dimensionSet dimPressure
tmp< fvMatrix< Type > > operator==(const fvMatrix< Type > &, const fvMatrix< Type > &)
const HashTable< dimensionSet > & dimensions()
Get the table of dimension sets.
Definition: dimensionSets.C:96
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
const dimensionSet dimDensity
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.