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 
31 #include "fvPatchFieldMapper.H"
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 {
75  patchType() = dict.getOrDefault<word>("patchType", word::null);
76 
77  if (dict.found("value"))
78  {
80  (
81  vectorField("value", dict, p.size())
82  );
83  }
84  else
85  {
86  fvPatchField<vector>::operator=(patchInternalField());
87  }
88 
89  refValue() = *this;
90  refGrad() = Zero;
91  valueFraction() = 1;
92 
93  if (pInf_ < SMALL)
94  {
96  << " unphysical pInf specified (pInf <= 0.0)"
97  << "\n on patch " << this->patch().name()
98  << " of field " << this->internalField().name()
99  << " in file " << this->internalField().objectPath()
100  << exit(FatalIOError);
101  }
102 }
103 
104 
107 (
109  const fvPatch& p,
111  const fvPatchFieldMapper& mapper
112 )
113 :
114  mixedFvPatchVectorField(ptf, p, iF, mapper),
115  TName_(ptf.TName_),
116  pName_(ptf.pName_),
117  psiName_(ptf.psiName_),
118  UInf_(ptf.UInf_),
119  pInf_(ptf.pInf_),
120  TInf_(ptf.TInf_),
121  gamma_(ptf.gamma_)
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  << "\n on patch " << this->patch().name()
188  << " of field " << this->internalField().name()
189  << " in file " << this->internalField().objectPath()
190  << exit(FatalError);
191  }
192 
193  vectorField& Up = refValue();
194  valueFraction() = 1;
195 
196  // get the near patch internal cell values
197  const vectorField U(patchInternalField());
198 
199 
200  // Find the component of U normal to the free-stream flow and in the
201  // plane of the free-stream flow and the patch normal
202 
203  // Direction of the free-stream flow
204  vector UInfHat = UInf_/mag(UInf_);
205 
206  // Normal to the plane defined by the free-stream and the patch normal
207  tmp<vectorField> nnInfHat = UInfHat ^ patch().nf();
208 
209  // Normal to the free-stream in the plane defined by the free-stream
210  // and the patch normal
211  const vectorField nHatInf(nnInfHat ^ UInfHat);
212 
213  // Component of U normal to the free-stream in the plane defined by the
214  // free-stream and the patch normal
215  const vectorField Un(nHatInf*(nHatInf & U));
216 
217  // The tangential component is
218  const vectorField Ut(U - Un);
219 
220  // Calculate the Prandtl-Meyer function of the free-stream
221  scalar nuMachInf =
222  sqrt((gamma_ + 1)/(gamma_ - 1))
223  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(MachInf) - 1)))
224  - atan(sqr(MachInf) - 1);
225 
226 
227  // Set the patch boundary condition based on the Mach number and direction
228  // of the flow dictated by the boundary/free-stream pressure difference
229 
230  forAll(Up, facei)
231  {
232  if (pp[facei] >= pInf_) // If outflow
233  {
234  // Assume supersonic outflow and calculate the boundary velocity
235  // according to ???
236 
237  scalar fpp =
238  sqrt(sqr(MachInf) - 1)
239  /(gamma_*sqr(MachInf))*mag(Ut[facei])*log(pp[facei]/pInf_);
240 
241  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
242 
243  // Calculate the Mach number of the boundary velocity
244  scalar Mach = mag(Up[facei])/sqrt(gamma_/ppsi[facei]);
245 
246  if (Mach <= 1) // If subsonic
247  {
248  // Zero-gradient subsonic outflow
249 
250  Up[facei] = U[facei];
251  valueFraction()[facei] = 0;
252  }
253  }
254  else // if inflow
255  {
256  // Calculate the Mach number of the boundary velocity
257  // from the boundary pressure assuming constant total pressure
258  // expansion from the free-stream
259  scalar Mach =
260  sqrt
261  (
262  (2/(gamma_ - 1))*(1 + ((gamma_ - 1)/2)*sqr(MachInf))
263  *pow(pp[facei]/pInf_, (1 - gamma_)/gamma_)
264  - 2/(gamma_ - 1)
265  );
266 
267  if (Mach > 1) // If supersonic
268  {
269  // Supersonic inflow is assumed to occur according to the
270  // Prandtl-Meyer expansion process
271 
272  scalar nuMachf =
273  sqrt((gamma_ + 1)/(gamma_ - 1))
274  *atan(sqrt((gamma_ - 1)/(gamma_ + 1)*(sqr(Mach) - 1)))
275  - atan(sqr(Mach) - 1);
276 
277  scalar fpp = (nuMachInf - nuMachf)*mag(Ut[facei]);
278 
279  Up[facei] = Ut[facei] + fpp*nHatInf[facei];
280  }
281  else // If subsonic
282  {
284  << "unphysical subsonic inflow has been generated"
285  << "\n on patch " << this->patch().name()
286  << " of field " << this->internalField().name()
287  << " in file "
288  << this->internalField().objectPath()
289  << exit(FatalError);
290  }
291  }
292  }
293 
294  mixedFvPatchVectorField::updateCoeffs();
295 }
296 
297 
299 {
301  os.writeEntryIfDifferent<word>("T", "T", TName_);
302  os.writeEntryIfDifferent<word>("p", "p", pName_);
303  os.writeEntryIfDifferent<word>("psi", "thermo:psi", psiName_);
304  os.writeEntry("UInf", UInf_);
305  os.writeEntry("pInf", pInf_);
306  os.writeEntry("TInf", TInf_);
307  os.writeEntry("gamma", gamma_);
308  writeEntry("value", os);
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 namespace Foam
315 {
317  (
320  );
321 }
322 
323 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: volSurfaceMapping.H:51
volFields.H
Foam::fvPatchField< vector >::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:384
Foam::Ostream::writeEntryIfDifferent
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:248
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::supersonicFreestreamFvPatchVectorField::write
virtual void write(Ostream &) const
Write.
Definition: supersonicFreestreamFvPatchVectorField.C:298
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::supersonicFreestreamFvPatchVectorField
This boundary condition provides a supersonic free-stream condition.
Definition: supersonicFreestreamFvPatchVectorField.H:128
Foam::FatalIOError
IOerror FatalIOError
fvPatchFieldMapper.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::vectorField
Field< vector > vectorField
Specialisation of Field<T> for vector.
Definition: primitiveFieldsFwd.H:54
supersonicFreestreamFvPatchVectorField.H
R
#define R(A, B, C, D, E, F, K, M)
Foam::Field< vector >
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:65
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::log
dimensionedScalar log(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:262
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::supersonicFreestreamFvPatchVectorField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: supersonicFreestreamFvPatchVectorField.C:162
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::supersonicFreestreamFvPatchVectorField::supersonicFreestreamFvPatchVectorField
supersonicFreestreamFvPatchVectorField(const fvPatch &, const DimensionedField< vector, volMesh > &)
Construct from patch and internal field.
Definition: supersonicFreestreamFvPatchVectorField.C:38
U
U
Definition: pEqn.H:72
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::Vector< scalar >
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::atan
dimensionedScalar atan(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:269
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:47
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::makePatchTypeField
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54