transonicEntrainmentPressureFvPatchScalarField.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) 2021-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 "surfaceFields.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
34 (
35  const fvPatch& p,
37  const dictionary& dict
38 )
39 :
40  mixedFvPatchScalarField(p, iF, dict, false),
41  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
42  psiName_(dict.lookupOrDefault<word>("psi", "psi")),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  gamma_(dict.lookup<scalar>("gamma")),
45  Mb_(dict.lookupOrDefault<scalar>("Mb", 0.5)),
46  p0_("p0", dict, p.size())
47 {
48  if (dict.found("value"))
49  {
51  (
52  scalarField("value", dict, p.size())
53  );
54  }
55  else
56  {
58  }
59 
60  refValue() = p0_;
61  refGrad() = Zero;
62  valueFraction() = 1;
63 }
64 
65 
68 (
70  const fvPatch& p,
72  const fvPatchFieldMapper& mapper
73 )
74 :
75  mixedFvPatchScalarField(psf, p, iF, mapper),
76  rhoName_(psf.rhoName_),
77  psiName_(psf.psiName_),
78  phiName_(psf.phiName_),
79  gamma_(psf.gamma_),
80  Mb_(psf.Mb_),
81  p0_(mapper(psf.p0_))
82 {}
83 
84 
87 (
90 )
91 :
92  mixedFvPatchScalarField(psf, iF),
93  rhoName_(psf.rhoName_),
94  psiName_(psf.psiName_),
95  phiName_(psf.phiName_),
96  gamma_(psf.gamma_),
97  Mb_(psf.Mb_),
98  p0_(psf.p0_)
99 {}
100 
101 
102 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
103 
105 (
106  const fvPatchScalarField& psf,
107  const fvPatchFieldMapper& mapper
108 )
109 {
110  mixedFvPatchScalarField::map(psf, mapper);
111 
113  refCast<const transonicEntrainmentPressureFvPatchScalarField>(psf);
114 
115  mapper(p0_, toppsf.p0_);
116 }
117 
118 
120 (
121  const fvPatchScalarField& psf
122 )
123 {
124  mixedFvPatchScalarField::reset(psf);
125 
127  refCast<const transonicEntrainmentPressureFvPatchScalarField>(psf);
128 
129  p0_.reset(toppsf.p0_);
130 }
131 
132 
134 {
135  if (updated())
136  {
137  return;
138  }
139 
140  const surfaceScalarField& phi =
141  db().lookupObject<surfaceScalarField>(phiName_);
142 
143  const fvsPatchField<scalar>& phip =
144  patch().patchField<surfaceScalarField, scalar>(phi);
145 
146  const fvPatchField<scalar>& psip =
147  patch().lookupPatchField<volScalarField, scalar>(psiName_);
148 
149  scalarField Unp(phip/patch().magSf());
150 
151  if (phi.dimensions() == dimMassFlux)
152  {
153  const fvPatchField<scalar>& rhop =
154  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
155 
156  Unp /= rhop;
157  }
158 
159  // Calculate the speed of sound and Mach number at the outlet patch
160  const scalarField c(sqrt(gamma_/psip));
161  const scalarField Ma(max(Unp/c, scalar(0)));
162 
163  const scalar gM1ByG = (gamma_ - 1)/gamma_;
164 
165  refValue() =
166  p0_
167  /pow
168  (
169  1 - (0.5*gM1ByG)*psip*negPart(Unp)*mag(Unp),
170  1/gM1ByG
171  );
172 
173  valueFraction() = 1 - min(max(Ma - Mb_, scalar(0))/(1 - Mb_), scalar(1));
174 
176 }
177 
178 
180 (
181  Ostream& os
182 ) const
183 {
185  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
186  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
187  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
188  writeEntry(os, "Mb", Mb_);
189  writeEntry(os, "gamma", gamma_);
190  writeEntry(os, "p0", p0_);
191  writeEntry(os, "value", *this);
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
196 
197 namespace Foam
198 {
200  (
203  );
204 }
205 
206 // ************************************************************************* //
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.
Generic GeometricField class.
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
friend Ostream & operator(Ostream &, const fvPatchField< Type > &)
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:204
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 an entrainment condition for pressure including support for superson...
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
transonicEntrainmentPressureFvPatchScalarField(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.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A class for handling words, derived from string.
Definition: word.H:62
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const dimensionSet dimMassFlux
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar negPart(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dictionary dict
volScalarField & p
Foam::surfaceFields.