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