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