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-2022 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 )
38 :
39  mixedFvPatchScalarField(p, iF),
40  rhoName_("rho"),
41  psiName_("thermo:psi"),
42  phiName_("phi"),
43  gamma_(0),
44  Mb_(0),
45  p0_(p.size(), 0)
46 {}
47 
48 
51 (
52  const fvPatch& p,
54  const dictionary& dict
55 )
56 :
57  mixedFvPatchScalarField(p, iF),
58  rhoName_(dict.lookupOrDefault<word>("rho", "rho")),
59  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
60  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
61  gamma_(dict.lookup<scalar>("gamma")),
62  Mb_(dict.lookupOrDefault<scalar>("Mb", 0.5)),
63  p0_("p0", dict, p.size())
64 {
65  if (dict.found("value"))
66  {
67  fvPatchScalarField::operator=
68  (
69  scalarField("value", dict, p.size())
70  );
71  }
72  else
73  {
74  fvPatchScalarField::operator=(p0_);
75  }
76 
77  refValue() = p0_;
78  refGrad() = Zero;
79  valueFraction() = 1;
80 }
81 
82 
85 (
87  const fvPatch& p,
89  const fvPatchFieldMapper& mapper
90 )
91 :
92  mixedFvPatchScalarField(psf, p, iF, mapper),
93  rhoName_(psf.rhoName_),
94  psiName_(psf.psiName_),
95  phiName_(psf.phiName_),
96  gamma_(psf.gamma_),
97  Mb_(psf.Mb_),
98  p0_(mapper(psf.p0_))
99 {}
100 
101 
104 (
107 )
108 :
109  mixedFvPatchScalarField(psf, iF),
110  rhoName_(psf.rhoName_),
111  psiName_(psf.psiName_),
112  phiName_(psf.phiName_),
113  gamma_(psf.gamma_),
114  Mb_(psf.Mb_),
115  p0_(psf.p0_)
116 {}
117 
118 
119 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
120 
122 (
123  const fvPatchFieldMapper& m
124 )
125 {
126  mixedFvPatchScalarField::autoMap(m);
127  m(p0_, p0_);
128 }
129 
130 
132 (
133  const fvPatchScalarField& psf,
134  const labelList& addr
135 )
136 {
137  mixedFvPatchScalarField::rmap(psf, addr);
138 
140  refCast<const transonicEntrainmentPressureFvPatchScalarField>(psf);
141 
142  p0_.rmap(toppsf.p0_, addr);
143 }
144 
145 
147 (
148  const fvPatchScalarField& psf
149 )
150 {
151  mixedFvPatchScalarField::reset(psf);
152 
154  refCast<const transonicEntrainmentPressureFvPatchScalarField>(psf);
155 
156  p0_.reset(toppsf.p0_);
157 }
158 
159 
161 {
162  if (updated())
163  {
164  return;
165  }
166 
167  const surfaceScalarField& phi =
168  db().lookupObject<surfaceScalarField>(phiName_);
169 
170  const fvsPatchField<scalar>& phip =
171  patch().patchField<surfaceScalarField, scalar>(phi);
172 
173  const fvPatchField<scalar>& psip =
174  patch().lookupPatchField<volScalarField, scalar>(psiName_);
175 
176  scalarField Unp(phip/patch().magSf());
177 
178  if (phi.dimensions() == dimMassFlux)
179  {
180  const fvPatchField<scalar>& rhop =
181  patch().lookupPatchField<volScalarField, scalar>(rhoName_);
182 
183  Unp /= rhop;
184  }
185 
186  // Calculate the speed of sound and Mach number at the outlet patch
187  const scalarField c(sqrt(gamma_/psip));
188  const scalarField Ma(max(Unp/c, scalar(0)));
189 
190  const scalar gM1ByG = (gamma_ - 1)/gamma_;
191 
192  refValue() =
193  p0_
194  /pow
195  (
196  1 - (0.5*gM1ByG)*psip*negPart(Unp)*mag(Unp),
197  1/gM1ByG
198  );
199 
200  valueFraction() = 1 - min(max(Ma - Mb_, scalar(0))/(1 - Mb_), scalar(1));
201 
203 }
204 
205 
207 (
208  Ostream& os
209 ) const
210 {
212  writeEntryIfDifferent<word>(os, "rho", "rho", rhoName_);
213  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
214  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
215  writeEntry(os, "Mb", Mb_);
216  writeEntry(os, "gamma", gamma_);
217  writeEntry(os, "p0", p0_);
218  writeEntry(os, "value", *this);
219 }
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 namespace Foam
225 {
227  (
230  );
231 }
232 
233 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:243
const dimensionedScalar c
Speed of light in a vacuum.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
This boundary condition provides an entrainment condition for pressure including support for superson...
const dimensionSet & dimensions() const
Return dimensions.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:157
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
volScalarField scalarField(fieldObject, mesh)
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimMassFlux
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: fvPatchField.C:216
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
transonicEntrainmentPressureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
Namespace for OpenFOAM.
dimensionedScalar negPart(const dimensionedScalar &ds)
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864