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 | Copyright (C) 2011-2016 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_(readScalar(dict.lookup("pInf"))),
69  TInf_(readScalar(dict.lookup("TInf"))),
70  gamma_(readScalar(dict.lookup("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 (
126 )
127 :
128  mixedFvPatchVectorField(sfspvf),
129  TName_(sfspvf.TName_),
130  pName_(sfspvf.pName_),
131  psiName_(sfspvf.psiName_),
132  UInf_(sfspvf.UInf_),
133  pInf_(sfspvf.pInf_),
134  TInf_(sfspvf.TInf_),
135  gamma_(sfspvf.gamma_)
136 {}
137 
138 
141 (
144 )
145 :
146  mixedFvPatchVectorField(sfspvf, iF),
147  TName_(sfspvf.TName_),
148  pName_(sfspvf.pName_),
149  psiName_(sfspvf.psiName_),
150  UInf_(sfspvf.UInf_),
151  pInf_(sfspvf.pInf_),
152  TInf_(sfspvf.TInf_),
153  gamma_(sfspvf.gamma_)
154 {}
155 
156 
157 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
158 
160 {
161  if (!size() || updated())
162  {
163  return;
164  }
165 
166  const fvPatchField<scalar>& pT =
167  patch().lookupPatchField<volScalarField, scalar>(TName_);
168 
169  const fvPatchField<scalar>& pp =
170  patch().lookupPatchField<volScalarField, scalar>(pName_);
171 
172  const fvPatchField<scalar>& ppsi =
173  patch().lookupPatchField<volScalarField, scalar>(psiName_);
174 
175  // Need R of the free-stream flow. Assume R is independent of location
176  // along patch so use face 0
177  scalar R = 1.0/(ppsi[0]*pT[0]);
178 
179  scalar MachInf = mag(UInf_)/sqrt(gamma_*R*TInf_);
180 
181  if (MachInf < 1.0)
182  {
184  << "\n on patch " << this->patch().name()
185  << " of field " << this->internalField().name()
186  << " in file " << this->internalField().objectPath()
187  << exit(FatalError);
188  }
189 
190  vectorField& Up = refValue();
191  valueFraction() = 1;
192 
193  // get the near patch internal cell values
194  const vectorField U(patchInternalField());
195 
196 
197  // Find the component of U normal to the free-stream flow and in the
198  // plane of the free-stream flow and the patch normal
199 
200  // Direction of the free-stream flow
201  vector UInfHat = UInf_/mag(UInf_);
202 
203  // Normal to the plane defined by the free-stream and the patch normal
204  tmp<vectorField> nnInfHat = UInfHat ^ patch().nf();
205 
206  // Normal to the free-stream in the plane defined by the free-stream
207  // and the patch normal
208  const vectorField nHatInf(nnInfHat ^ UInfHat);
209 
210  // Component of U normal to the free-stream in the plane defined by the
211  // free-stream and the patch normal
212  const vectorField Un(nHatInf*(nHatInf & U));
213 
214  // The tangential component is
215  const vectorField Ut(U - Un);
216 
217  // Calculate the Prandtl-Meyer function of the free-stream
218  scalar nuMachInf =
219  sqrt((gamma_ + 1)/(gamma_ - 1))
220  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(MachInf) - 1)))
221  - atan(sqr(MachInf) - 1);
222 
223 
224  // Set the patch boundary condition based on the Mach number and direction
225  // of the flow dictated by the boundary/free-stream pressure difference
226 
227  forAll(Up, facei)
228  {
229  if (pp[facei] >= pInf_) // If outflow
230  {
231  // Assume supersonic outflow and calculate the boundary velocity
232  // according to ???
233 
234  scalar fpp =
235  sqrt(sqr(MachInf) - 1)
236  /(gamma_*sqr(MachInf))*mag(Ut[facei])*log(pp[facei]/pInf_);
237 
238  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
239 
240  // Calculate the Mach number of the boundary velocity
241  scalar Mach = mag(Up[facei])/sqrt(gamma_/ppsi[facei]);
242 
243  if (Mach <= 1) // If subsonic
244  {
245  // Zero-gradient subsonic outflow
246 
247  Up[facei] = U[facei];
248  valueFraction()[facei] = 0;
249  }
250  }
251  else // if inflow
252  {
253  // Calculate the Mach number of the boundary velocity
254  // from the boundary pressure assuming constant total pressure
255  // expansion from the free-stream
256  scalar Mach =
257  sqrt
258  (
259  (2/(gamma_ - 1))*(1 + ((gamma_ - 1)/2)*sqr(MachInf))
260  *pow(pp[facei]/pInf_, (1 - gamma_)/gamma_)
261  - 2/(gamma_ - 1)
262  );
263 
264  if (Mach > 1) // If supersonic
265  {
266  // Supersonic inflow is assumed to occur according to the
267  // Prandtl-Meyer expansion process
268 
269  scalar nuMachf =
270  sqrt((gamma_ + 1)/(gamma_ - 1))
271  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(Mach) - 1)))
272  - atan(sqr(Mach) - 1);
273 
274  scalar fpp = (nuMachInf - nuMachf)*mag(Ut[facei]);
275 
276  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
277  }
278  else // If subsonic
279  {
281  << "unphysical subsonic inflow has been generated"
282  << "\n on patch " << this->patch().name()
283  << " of field " << this->internalField().name()
284  << " in file "
285  << this->internalField().objectPath()
286  << exit(FatalError);
287  }
288  }
289  }
290 
291  mixedFvPatchVectorField::updateCoeffs();
292 }
293 
294 
296 {
298  writeEntryIfDifferent<word>(os, "T", "T", TName_);
299  writeEntryIfDifferent<word>(os, "p", "p", pName_);
300  writeEntryIfDifferent<word>(os, "psi", "thermo:psi", psiName_);
301  os.writeKeyword("UInf") << UInf_ << token::END_STATEMENT << nl;
302  os.writeKeyword("pInf") << pInf_ << token::END_STATEMENT << nl;
303  os.writeKeyword("TInf") << TInf_ << token::END_STATEMENT << nl;
304  os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
305  writeEntry("value", os);
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 namespace Foam
312 {
314  (
317  );
318 }
319 
320 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
U
Definition: pEqn.H:83
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:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
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:362
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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:91
virtual label size() const
Return size.
Definition: fvPatch.H:161
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
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:331
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
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:576
IOerror FatalIOError