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-2023 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 "fvPatchFieldMapper.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") : 1),
48  p0_(Function1<scalar>::New("p0", dict))
49 {
50  if (dict.found("value"))
51  {
53  (
54  scalarField("value", dict, p.size())
55  );
56  }
57  else
58  {
59  const scalar t = this->db().time().userTimeValue();
60  fvPatchScalarField::operator==(p0_->value(t));
61  }
62 }
63 
64 
67 (
69  const fvPatch& p,
71  const fvPatchFieldMapper& mapper
72 )
73 :
74  fixedValueFvPatchScalarField(ptf, p, iF, mapper, false), // Don't map
75  UName_(ptf.UName_),
76  phiName_(ptf.phiName_),
77  rhoName_(ptf.rhoName_),
78  psiName_(ptf.psiName_),
79  gamma_(ptf.gamma_),
80  p0_(ptf.p0_, false)
81 {
82  // Set the patch pressure to the current total pressure
83  // This is not ideal but avoids problems with the creation of patch faces
84  const scalar t = this->db().time().userTimeValue();
85  fvPatchScalarField::operator==(p0_->value(t));
86 }
87 
88 
91 (
94 )
95 :
96  fixedValueFvPatchScalarField(ptf, iF),
97  UName_(ptf.UName_),
98  phiName_(ptf.phiName_),
99  rhoName_(ptf.rhoName_),
100  psiName_(ptf.psiName_),
101  gamma_(ptf.gamma_),
102  p0_(ptf.p0_, false)
103 {}
104 
105 
106 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
107 
109 (
110  const vectorField& Up
111 )
112 {
113  if (updated())
114  {
115  return;
116  }
117 
118  scalar p0 = p0_->value(this->db().time().userTimeValue());
119 
120  const fvsPatchField<scalar>& phip =
121  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
122 
123  if (internalField().dimensions() == dimPressure)
124  {
125  if (psiName_ == "none")
126  {
127  // Variable density and low-speed compressible flow
128 
129  const fvPatchField<scalar>& rho =
130  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
131 
132  operator==(p0 - 0.5*rho*neg(phip)*magSqr(Up));
133  }
134  else
135  {
136  // High-speed compressible flow
137 
138  const fvPatchField<scalar>& psip =
139  patch().lookupPatchField<volScalarField, scalar>(psiName_);
140 
141  if (gamma_ > 1)
142  {
143  scalar gM1ByG = (gamma_ - 1)/gamma_;
144 
145  operator==
146  (
147  p0
148  /pow
149  (
150  (1.0 + 0.5*psip*gM1ByG*neg(phip)*magSqr(Up)),
151  1.0/gM1ByG
152  )
153  );
154  }
155  else
156  {
157  operator==(p0/(1.0 + 0.5*psip*neg(phip)*magSqr(Up)));
158  }
159  }
160 
161  }
162  else if (internalField().dimensions() == dimPressure/dimDensity)
163  {
164  // Incompressible flow
165  operator==(p0 - 0.5*neg(phip)*magSqr(Up));
166  }
167  else
168  {
170  << " Incorrect pressure dimensions " << internalField().dimensions()
171  << nl
172  << " Should be " << dimPressure
173  << " for compressible/variable density flow" << nl
174  << " or " << dimPressure/dimDensity
175  << " for incompressible flow," << nl
176  << " on patch " << this->patch().name()
177  << " of field " << this->internalField().name()
178  << " in file " << this->internalField().objectPath()
179  << exit(FatalError);
180  }
181 
182  fixedValueFvPatchScalarField::updateCoeffs();
183 }
184 
185 
187 {
188  updateCoeffs(patch().lookupPatchField<volVectorField, vector>(UName_));
189 }
190 
191 
193 {
195  writeEntryIfDifferent<word>(os, "U", "U", UName_);
196  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
197  writeEntry(os, "rho", rhoName_);
198  writeEntry(os, "psi", psiName_);
199  writeEntry(os, "gamma", gamma_);
200  writeEntry(os, p0_());
201  writeEntry(os, "value", *this);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
207 namespace Foam
208 {
210  (
213  );
214 }
215 
216 // ************************************************************************* //
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...
Run-time selectable general function of one variable.
Definition: Function1.H:64
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:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator==(const fvPatchField< Type > &)
Definition: fvPatchField.C:417
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:82
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:306
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 > &)
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:61
dimensionedScalar neg(const dimensionedScalar &ds)
error FatalError
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
static const char nl
Definition: Ostream.H:260
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.