supersonicFreestreamFvPatchVectorField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2017-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
32 #include "volFields.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
38 (
39  const fvPatch& p,
41 )
42 :
43  mixedFvPatchVectorField(p, iF),
44  TName_("T"),
45  pName_("p"),
46  psiName_("thermo:psi"),
47  UInf_(Zero),
48  pInf_(0),
49  TInf_(0),
50  gamma_(0)
51 {
52  refValue() = patchInternalField();
53  refGrad() = Zero;
54  valueFraction() = 1;
55 }
56 
57 
60 (
61  const fvPatch& p,
63  const dictionary& dict
64 )
65 :
66  mixedFvPatchVectorField(p, iF),
67  TName_(dict.getOrDefault<word>("T", "T")),
68  pName_(dict.getOrDefault<word>("p", "p")),
69  psiName_(dict.getOrDefault<word>("psi", "thermo:psi")),
70  UInf_(dict.lookup("UInf")),
71  pInf_(dict.get<scalar>("pInf")),
72  TInf_(dict.get<scalar>("TInf")),
73  gamma_(dict.get<scalar>("gamma"))
74 {
76 
77  if (!this->readValueEntry(dict))
78  {
79  this->extrapolateInternal();
80  }
81 
82  refValue() = *this;
83  refGrad() = Zero;
84  valueFraction() = 1;
85 
86  if (pInf_ < SMALL)
87  {
89  << " unphysical pInf specified (pInf <= 0.0)"
90  << "\n on patch " << this->patch().name()
91  << " of field " << this->internalField().name()
92  << " in file " << this->internalField().objectPath()
94  }
95 }
96 
97 
100 (
101  const supersonicFreestreamFvPatchVectorField& ptf,
102  const fvPatch& p,
103  const DimensionedField<vector, volMesh>& iF,
104  const fvPatchFieldMapper& mapper
105 )
106 :
107  mixedFvPatchVectorField(ptf, p, iF, mapper),
108  TName_(ptf.TName_),
109  pName_(ptf.pName_),
110  psiName_(ptf.psiName_),
111  UInf_(ptf.UInf_),
112  pInf_(ptf.pInf_),
113  TInf_(ptf.TInf_),
114  gamma_(ptf.gamma_)
115 {}
116 
117 
120 (
122 )
123 :
124  mixedFvPatchVectorField(sfspvf),
125  TName_(sfspvf.TName_),
126  pName_(sfspvf.pName_),
127  psiName_(sfspvf.psiName_),
128  UInf_(sfspvf.UInf_),
129  pInf_(sfspvf.pInf_),
130  TInf_(sfspvf.TInf_),
131  gamma_(sfspvf.gamma_)
132 {}
133 
134 
137 (
140 )
141 :
142  mixedFvPatchVectorField(sfspvf, iF),
143  TName_(sfspvf.TName_),
144  pName_(sfspvf.pName_),
145  psiName_(sfspvf.psiName_),
146  UInf_(sfspvf.UInf_),
147  pInf_(sfspvf.pInf_),
148  TInf_(sfspvf.TInf_),
149  gamma_(sfspvf.gamma_)
150 {}
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
156 {
157  if (!size() || updated())
158  {
159  return;
160  }
161 
162  const auto& pT = patch().lookupPatchField<volScalarField>(TName_);
163 
164  const auto& pp = patch().lookupPatchField<volScalarField>(pName_);
165 
166  const auto& ppsi = patch().lookupPatchField<volScalarField>(psiName_);
167 
168  // Need R of the free-stream flow. Assume R is independent of location
169  // along patch so use face 0
170  scalar R = 1.0/(ppsi[0]*pT[0]);
171 
172  scalar MachInf = mag(UInf_)/sqrt(gamma_*R*TInf_);
173 
174  if (MachInf < 1.0)
175  {
177  << "\n on patch " << this->patch().name()
178  << " of field " << this->internalField().name()
179  << " in file " << this->internalField().objectPath()
180  << exit(FatalError);
181  }
182 
183  vectorField& Up = refValue();
184  valueFraction() = 1;
185 
186  // get the near patch internal cell values
187  const vectorField U(patchInternalField());
188 
189 
190  // Find the component of U normal to the free-stream flow and in the
191  // plane of the free-stream flow and the patch normal
192 
193  // Direction of the free-stream flow
194  vector UInfHat = UInf_/mag(UInf_);
195 
196  // Normal to the plane defined by the free-stream and the patch normal
197  tmp<vectorField> nnInfHat = UInfHat ^ patch().nf();
198 
199  // Normal to the free-stream in the plane defined by the free-stream
200  // and the patch normal
201  const vectorField nHatInf(nnInfHat ^ UInfHat);
202 
203  // Component of U normal to the free-stream in the plane defined by the
204  // free-stream and the patch normal
205  const vectorField Un(nHatInf*(nHatInf & U));
206 
207  // The tangential component is
208  const vectorField Ut(U - Un);
209 
210  // Calculate the Prandtl-Meyer function of the free-stream
211  scalar nuMachInf =
212  sqrt((gamma_ + 1)/(gamma_ - 1))
213  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(MachInf) - 1)))
214  - atan(sqr(MachInf) - 1);
215 
216 
217  // Set the patch boundary condition based on the Mach number and direction
218  // of the flow dictated by the boundary/free-stream pressure difference
219 
220  forAll(Up, facei)
221  {
222  if (pp[facei] >= pInf_) // If outflow
223  {
224  // Assume supersonic outflow and calculate the boundary velocity
225  // according to ???
226 
227  scalar fpp =
228  sqrt(sqr(MachInf) - 1)
229  /(gamma_*sqr(MachInf))*mag(Ut[facei])*log(pp[facei]/pInf_);
230 
231  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
232 
233  // Calculate the Mach number of the boundary velocity
234  scalar Mach = mag(Up[facei])/sqrt(gamma_/ppsi[facei]);
235 
236  if (Mach <= 1) // If subsonic
237  {
238  // Zero-gradient subsonic outflow
239 
240  Up[facei] = U[facei];
241  valueFraction()[facei] = 0;
242  }
243  }
244  else // if inflow
245  {
246  // Calculate the Mach number of the boundary velocity
247  // from the boundary pressure assuming constant total pressure
248  // expansion from the free-stream
249  scalar Mach =
250  sqrt
251  (
252  (2/(gamma_ - 1))*(1 + ((gamma_ - 1)/2)*sqr(MachInf))
253  *pow(pp[facei]/pInf_, (1 - gamma_)/gamma_)
254  - 2/(gamma_ - 1)
255  );
256 
257  if (Mach > 1) // If supersonic
258  {
259  // Supersonic inflow is assumed to occur according to the
260  // Prandtl-Meyer expansion process
261 
262  scalar nuMachf =
263  sqrt((gamma_ + 1)/(gamma_ - 1))
264  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(Mach) - 1)))
265  - atan(sqr(Mach) - 1);
266 
267  scalar fpp = (nuMachInf - nuMachf)*mag(Ut[facei]);
268 
269  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
270  }
271  else // If subsonic
272  {
274  << "unphysical subsonic inflow has been generated"
275  << "\n on patch " << this->patch().name()
276  << " of field " << this->internalField().name()
277  << " in file "
278  << this->internalField().objectPath()
279  << exit(FatalError);
280  }
281  }
282  }
283 
284  mixedFvPatchVectorField::updateCoeffs();
285 }
286 
287 
289 {
291  os.writeEntryIfDifferent<word>("T", "T", TName_);
292  os.writeEntryIfDifferent<word>("p", "p", pName_);
293  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
294  os.writeEntry("UInf", UInf_);
295  os.writeEntry("pInf", pInf_);
296  os.writeEntry("TInf", TInf_);
297  os.writeEntry("gamma", gamma_);
299 }
300 
301 
302 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
303 
304 namespace Foam
305 {
307  (
309  supersonicFreestreamFvPatchVectorField
310  );
311 }
312 
313 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
dictionary dict
fvPatchField< vector > fvPatchVectorField
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
Definition: fvPatchField.H:375
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:372
virtual void readDict(const dictionary &dict)
Read dictionary entries.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Vector< scalar > vector
Definition: vector.H:57
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
OBJstream os(runTime.globalPath()/outputName)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
This boundary condition provides a supersonic free-stream condition.
#define R(A, B, C, D, E, F, K, M)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
dimensionedScalar atan(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
Field< vector > vectorField
Specialisation of Field<T> for vector.
volScalarField & p
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
supersonicFreestreamFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127