inletOutletTotalTemperatureFvPatchScalarField.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-2023 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 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
36 (
37  const fvPatch& p,
39  const dictionary& dict
40 )
41 :
42  inletOutletFvPatchScalarField(p, iF),
43  UName_(dict.lookupOrDefault<word>("U", "U")),
44  psiName_(dict.lookupOrDefault<word>("psi", "psi")),
45  gamma_(dict.lookup<scalar>("gamma")),
46  T0_("T0", dict, p.size())
47 {
48  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
49 
50  this->refValue() = Zero;
51  if (dict.found("value"))
52  {
54  (
55  scalarField("value", dict, p.size())
56  );
57  }
58  else
59  {
61  }
62 
63  this->refGrad() = Zero;
64  this->valueFraction() = 0.0;
65 }
66 
67 
70 (
72  const fvPatch& p,
74  const fvPatchFieldMapper& mapper
75 )
76 :
77  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
78  UName_(ptf.UName_),
79  psiName_(ptf.psiName_),
80  gamma_(ptf.gamma_),
81  T0_(mapper(ptf.T0_))
82 {}
83 
84 
87 (
90 )
91 :
92  inletOutletFvPatchScalarField(tppsf, iF),
93  UName_(tppsf.UName_),
94  psiName_(tppsf.psiName_),
95  gamma_(tppsf.gamma_),
96  T0_(tppsf.T0_)
97 {}
98 
99 
100 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
101 
103 (
104  const fvPatchScalarField& ptf,
105  const fvPatchFieldMapper& mapper
106 )
107 {
108  inletOutletFvPatchScalarField::map(ptf, mapper);
109 
111  refCast<const inletOutletTotalTemperatureFvPatchScalarField>(ptf);
112 
113  mapper(T0_, tiptf.T0_);
114 }
115 
116 
118 (
119  const fvPatchScalarField& ptf
120 )
121 {
122  inletOutletFvPatchScalarField::reset(ptf);
123 
125  refCast<const inletOutletTotalTemperatureFvPatchScalarField>(ptf);
126 
127  T0_.reset(tiptf.T0_);
128 }
129 
130 
132 {
133  if (updated())
134  {
135  return;
136  }
137 
138  const fvPatchVectorField& Up =
139  patch().lookupPatchField<volVectorField, vector>(UName_);
140 
141  const fvsPatchField<scalar>& phip =
142  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
143 
144  const fvPatchField<scalar>& psip =
145  patch().lookupPatchField<volScalarField, scalar>(psiName_);
146 
147  scalar gM1ByG = (gamma_ - 1.0)/gamma_;
148 
149  this->refValue() =
150  T0_/(1.0 + 0.5*psip*gM1ByG*neg(phip)*magSqr(Up));
151  this->valueFraction() = neg(phip);
152 
153  inletOutletFvPatchScalarField::updateCoeffs();
154 }
155 
156 
158 const
159 {
161  writeEntryIfDifferent<word>(os, "U", "U", UName_);
162  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
163  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
164  writeEntry(os, "gamma", gamma_);
165  writeEntry(os, "T0", T0_);
166  writeEntry(os, "value", *this);
167 }
168 
169 
170 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
171 
172 namespace Foam
173 {
175  (
178  );
179 }
180 
181 // ************************************************************************* //
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...
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:160
Foam::fvPatchFieldMapper.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:87
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:231
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:251
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
This boundary condition provides an outflow condition for total temperature for use with supersonic c...
inletOutletTotalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
virtual void reset(const fvPatchScalarField &)
Reset the fvPatchField to the given fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
virtual void map(const fvPatchScalarField &, const fvPatchFieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
SurfaceField< scalar > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:61
dimensionedScalar neg(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.