dynamicPressureFvPatchScalarField.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 
27 #include "fvPatchFieldMapper.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
37 }
38 
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
43 (
44  const scalarField& p0p,
45  const scalarField& K0mKp
46 )
47 {
48  if (updated())
49  {
50  return;
51  }
52 
53  if (internalField().dimensions() == dimPressure)
54  {
55  if (psiName_ == "none")
56  {
57  // Variable density and low-speed compressible flow
58 
59  const fvPatchField<scalar>& rho =
60  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
61 
62  operator==(p0p + rho*K0mKp);
63  }
64  else
65  {
66  // High-speed compressible flow
67 
68  const fvPatchField<scalar>& psip =
69  patch().lookupPatchField<volScalarField, scalar>(psiName_);
70 
71  if (gamma_ > 1)
72  {
73  const scalar gM1ByG = (gamma_ - 1)/gamma_;
74 
75  operator==
76  (
77  p0p/pow(scalar(1) - psip*gM1ByG*K0mKp, 1/gM1ByG)
78  );
79  }
80  else
81  {
82  operator==(p0p/(scalar(1) - psip*K0mKp));
83  }
84  }
85  }
86  else if (internalField().dimensions() == dimPressure/dimDensity)
87  {
88  // Incompressible flow
89 
90  operator==(p0p + K0mKp);
91  }
92  else
93  {
95  << " Incorrect pressure dimensions " << internalField().dimensions()
96  << nl
97  << " Should be " << dimPressure
98  << " for compressible/variable density flow" << nl
99  << " or " << dimPressure/dimDensity
100  << " for incompressible flow," << nl
101  << " on patch " << this->patch().name()
102  << " of field " << this->internalField().name()
103  << " in file " << this->internalField().objectPath()
104  << exit(FatalError);
105  }
106 
107  fixedValueFvPatchScalarField::updateCoeffs();
108 }
109 
110 
111 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
112 
114 (
115  const fvPatch& p,
117  const dictionary& dict
118 )
119 :
120  fixedValueFvPatchScalarField(p, iF, dict, false),
121  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
122  psiName_(dict.lookupOrDefault<word>("psi", "none")),
123  gamma_(dict.lookupOrDefault<scalar>("gamma", 1)),
124  p0_("p0", dict, p.size())
125 {
126  if (dict.found("value"))
127  {
129  (
130  scalarField("value", dict, p.size())
131  );
132  }
133  else
134  {
136  }
137 }
138 
139 
141 (
143  const fvPatch& p,
145  const fvPatchFieldMapper& mapper
146 )
147 :
148  fixedValueFvPatchScalarField(psf, p, iF, mapper),
149  rhoName_(psf.rhoName_),
150  psiName_(psf.psiName_),
151  gamma_(psf.gamma_),
152  p0_(mapper(psf.p0_))
153 {}
154 
155 
157 (
160 )
161 :
162  fixedValueFvPatchScalarField(tppsf, iF),
163  rhoName_(tppsf.rhoName_),
164  psiName_(tppsf.psiName_),
165  gamma_(tppsf.gamma_),
166  p0_(tppsf.p0_)
167 {}
168 
169 
170 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
171 
173 (
174  const fvPatchScalarField& psf,
175  const fvPatchFieldMapper& mapper
176 )
177 {
178  fixedValueFvPatchScalarField::map(psf, mapper);
179 
181  refCast<const dynamicPressureFvPatchScalarField>(psf);
182 
183  mapper(p0_, dpsf.p0_);
184 }
185 
186 
188 (
189  const fvPatchScalarField& psf
190 )
191 {
192  fixedValueFvPatchScalarField::reset(psf);
193 
195  refCast<const dynamicPressureFvPatchScalarField>(psf);
196 
197  p0_.reset(dpsf.p0_);
198 }
199 
200 
202 {
204  writeEntry(os, "rho", rhoName_);
205  writeEntry(os, "psi", psiName_);
206  writeEntry(os, "gamma", gamma_);
207  writeEntry(os, "p0", p0_);
208  writeEntry(os, "value", *this);
209 }
210 
211 
212 // ************************************************************************* //
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:57
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
This boundary condition provides a dynamic pressure condition. It subtracts a kinetic energy term fro...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
void updateCoeffs(const scalarField &p0p, const scalarField &K0mKp)
Update the coefficients associated with the patch field.
dynamicPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
const word rhoName_
Name of the density field used to normalise the mass flux.
const word psiName_
Name of the compressibility field used to calculate the wave speed.
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 UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
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
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 > &)
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
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
static const char nl
Definition: Ostream.H:260
dictionary dict
volScalarField & p
Foam::surfaceFields.