supersonicFreestreamFvPatchVectorField.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-2020 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 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
35 (
36  const fvPatch& p,
38 )
39 :
40  mixedFvPatchVectorField(p, iF),
41  TName_("T"),
42  pName_("p"),
43  psiName_("thermo:psi"),
44  UInf_(Zero),
45  pInf_(0),
46  TInf_(0),
47  gamma_(0)
48 {
49  refValue() = patchInternalField();
50  refGrad() = Zero;
51  valueFraction() = 1;
52 }
53 
54 
57 (
58  const fvPatch& p,
60  const dictionary& dict
61 )
62 :
63  mixedFvPatchVectorField(p, iF),
64  TName_(dict.lookupOrDefault<word>("T", "T")),
65  pName_(dict.lookupOrDefault<word>("p", "p")),
66  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
67  UInf_(dict.lookup("UInf")),
68  pInf_(dict.lookup<scalar>("pInf")),
69  TInf_(dict.lookup<scalar>("TInf")),
70  gamma_(dict.lookup<scalar>("gamma"))
71 {
72  if (dict.found("value"))
73  {
75  (
76  vectorField("value", dict, p.size())
77  );
78  }
79  else
80  {
81  fvPatchField<vector>::operator=(patchInternalField());
82  }
83 
84  refValue() = *this;
85  refGrad() = Zero;
86  valueFraction() = 1;
87 
88  if (pInf_ < small)
89  {
91  (
92  dict
93  ) << " unphysical pInf specified (pInf <= 0.0)"
94  << "\n on patch " << this->patch().name()
95  << " of field " << this->internalField().name()
96  << " in file " << this->internalField().objectPath()
97  << exit(FatalIOError);
98  }
99 }
100 
101 
104 (
106  const fvPatch& p,
108  const fvPatchFieldMapper& mapper
109 )
110 :
111  mixedFvPatchVectorField(ptf, p, iF, mapper),
112  TName_(ptf.TName_),
113  pName_(ptf.pName_),
114  psiName_(ptf.psiName_),
115  UInf_(ptf.UInf_),
116  pInf_(ptf.pInf_),
117  TInf_(ptf.TInf_),
118  gamma_(ptf.gamma_)
119 {}
120 
121 
124 (
127 )
128 :
129  mixedFvPatchVectorField(sfspvf, iF),
130  TName_(sfspvf.TName_),
131  pName_(sfspvf.pName_),
132  psiName_(sfspvf.psiName_),
133  UInf_(sfspvf.UInf_),
134  pInf_(sfspvf.pInf_),
135  TInf_(sfspvf.TInf_),
136  gamma_(sfspvf.gamma_)
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
143 {
144  if (!size() || updated())
145  {
146  return;
147  }
148 
149  const fvPatchField<scalar>& pT =
150  patch().lookupPatchField<volScalarField, scalar>(TName_);
151 
152  const fvPatchField<scalar>& pp =
153  patch().lookupPatchField<volScalarField, scalar>(pName_);
154 
155  const fvPatchField<scalar>& ppsi =
156  patch().lookupPatchField<volScalarField, scalar>(psiName_);
157 
158  // Need R of the free-stream flow. Assume R is independent of location
159  // along patch so use face 0
160  scalar R = 1.0/(ppsi[0]*pT[0]);
161 
162  scalar MachInf = mag(UInf_)/sqrt(gamma_*R*TInf_);
163 
164  if (MachInf < 1.0)
165  {
167  << "\n on patch " << this->patch().name()
168  << " of field " << this->internalField().name()
169  << " in file " << this->internalField().objectPath()
170  << exit(FatalError);
171  }
172 
173  vectorField& Up = refValue();
174  valueFraction() = 1;
175 
176  // get the near patch internal cell values
177  const vectorField U(patchInternalField());
178 
179 
180  // Find the component of U normal to the free-stream flow and in the
181  // plane of the free-stream flow and the patch normal
182 
183  // Direction of the free-stream flow
184  vector UInfHat = UInf_/mag(UInf_);
185 
186  // Normal to the plane defined by the free-stream and the patch normal
187  tmp<vectorField> nnInfHat = UInfHat ^ patch().nf();
188 
189  // Normal to the free-stream in the plane defined by the free-stream
190  // and the patch normal
191  const vectorField nHatInf(nnInfHat ^ UInfHat);
192 
193  // Component of U normal to the free-stream in the plane defined by the
194  // free-stream and the patch normal
195  const vectorField Un(nHatInf*(nHatInf & U));
196 
197  // The tangential component is
198  const vectorField Ut(U - Un);
199 
200  // Calculate the Prandtl-Meyer function of the free-stream
201  scalar nuMachInf =
202  sqrt((gamma_ + 1)/(gamma_ - 1))
203  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(MachInf) - 1)))
204  - atan(sqr(MachInf) - 1);
205 
206 
207  // Set the patch boundary condition based on the Mach number and direction
208  // of the flow dictated by the boundary/free-stream pressure difference
209 
210  forAll(Up, facei)
211  {
212  if (pp[facei] >= pInf_) // If outflow
213  {
214  // Assume supersonic outflow and calculate the boundary velocity
215  // according to ???
216 
217  scalar fpp =
218  sqrt(sqr(MachInf) - 1)
219  /(gamma_*sqr(MachInf))*mag(Ut[facei])*log(pp[facei]/pInf_);
220 
221  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
222 
223  // Calculate the Mach number of the boundary velocity
224  scalar Mach = mag(Up[facei])/sqrt(gamma_/ppsi[facei]);
225 
226  if (Mach <= 1) // If subsonic
227  {
228  // Zero-gradient subsonic outflow
229 
230  Up[facei] = U[facei];
231  valueFraction()[facei] = 0;
232  }
233  }
234  else // if inflow
235  {
236  // Calculate the Mach number of the boundary velocity
237  // from the boundary pressure assuming constant total pressure
238  // expansion from the free-stream
239  scalar Mach =
240  sqrt
241  (
242  (2/(gamma_ - 1))*(1 + ((gamma_ - 1)/2)*sqr(MachInf))
243  *pow(pp[facei]/pInf_, (1 - gamma_)/gamma_)
244  - 2/(gamma_ - 1)
245  );
246 
247  if (Mach > 1) // If supersonic
248  {
249  // Supersonic inflow is assumed to occur according to the
250  // Prandtl-Meyer expansion process
251 
252  scalar nuMachf =
253  sqrt((gamma_ + 1)/(gamma_ - 1))
254  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(Mach) - 1)))
255  - atan(sqr(Mach) - 1);
256 
257  scalar fpp = (nuMachInf - nuMachf)*mag(Ut[facei]);
258 
259  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
260  }
261  else // If subsonic
262  {
264  << "unphysical subsonic inflow has been generated"
265  << "\n on patch " << this->patch().name()
266  << " of field " << this->internalField().name()
267  << " in file "
268  << this->internalField().objectPath()
269  << exit(FatalError);
270  }
271  }
272  }
273 
274  mixedFvPatchVectorField::updateCoeffs();
275 }
276 
277 
279 {
281  writeEntryIfDifferent<word>(os, "T", "T", TName_);
282  writeEntryIfDifferent<word>(os, "p", "p", pName_);
283  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
284  writeEntry(os, "UInf", UInf_);
285  writeEntry(os, "pInf", pInf_);
286  writeEntry(os, "TInf", TInf_);
287  writeEntry(os, "gamma", gamma_);
288  writeEntry(os, "value", *this);
289 }
290 
291 
292 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
293 
294 namespace Foam
295 {
297  (
300  );
301 }
302 
303 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:62
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
volVectorField vectorField(fieldObject, mesh)
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:260
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
static const zero Zero
Definition: zero.H:97
virtual label size() const
Return size.
Definition: fvPatch.H:156
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
This boundary condition provides a supersonic free-stream condition.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define R(A, B, C, D, E, F, K, M)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:335
dimensionedScalar atan(const dimensionedScalar &ds)
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
supersonicFreestreamFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844
IOerror FatalIOError