inletOutletTotalTemperatureFvPatchScalarField.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-2017 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 )
40 :
41  inletOutletFvPatchScalarField(p, iF),
42  UName_("U"),
43  psiName_("psi"),
44  gamma_(0.0),
45  T0_(p.size(), 0.0)
46 {
47  this->refValue() = Zero;
48  this->refGrad() = Zero;
49  this->valueFraction() = 0.0;
50 }
51 
52 
55 (
57  const fvPatch& p,
59  const fvPatchFieldMapper& mapper
60 )
61 :
62  inletOutletFvPatchScalarField(ptf, p, iF, mapper),
63  UName_(ptf.UName_),
64  psiName_(ptf.psiName_),
65  gamma_(ptf.gamma_),
66  T0_(ptf.T0_, mapper)
67 {}
68 
69 
72 (
73  const fvPatch& p,
75  const dictionary& dict
76 )
77 :
78  inletOutletFvPatchScalarField(p, iF),
79  UName_(dict.lookupOrDefault<word>("U", "U")),
80  psiName_(dict.lookupOrDefault<word>("psi", "thermo:psi")),
81  gamma_(readScalar(dict.lookup("gamma"))),
82  T0_("T0", dict, p.size())
83 {
84  this->phiName_ = dict.lookupOrDefault<word>("phi", "phi");
85 
86  this->refValue() = Zero;
87  if (dict.found("value"))
88  {
90  (
91  scalarField("value", dict, p.size())
92  );
93  }
94  else
95  {
97  }
98 
99  this->refGrad() = Zero;
100  this->valueFraction() = 0.0;
101 }
102 
103 
106 (
108 )
109 :
110  inletOutletFvPatchScalarField(tppsf),
111  UName_(tppsf.UName_),
112  psiName_(tppsf.psiName_),
113  gamma_(tppsf.gamma_),
114  T0_(tppsf.T0_)
115 {}
116 
117 
120 (
123 )
124 :
125  inletOutletFvPatchScalarField(tppsf, iF),
126  UName_(tppsf.UName_),
127  psiName_(tppsf.psiName_),
128  gamma_(tppsf.gamma_),
129  T0_(tppsf.T0_)
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
136 (
137  const fvPatchFieldMapper& m
138 )
139 {
140  inletOutletFvPatchScalarField::autoMap(m);
141  T0_.autoMap(m);
142 }
143 
144 
146 (
147  const fvPatchScalarField& ptf,
148  const labelList& addr
149 )
150 {
151  inletOutletFvPatchScalarField::rmap(ptf, addr);
152 
154  refCast<const inletOutletTotalTemperatureFvPatchScalarField>(ptf);
155 
156  T0_.rmap(tiptf.T0_, addr);
157 }
158 
159 
161 {
162  if (updated())
163  {
164  return;
165  }
166 
167  const fvPatchVectorField& Up =
168  patch().lookupPatchField<volVectorField, vector>(UName_);
169 
170  const fvsPatchField<scalar>& phip =
171  patch().lookupPatchField<surfaceScalarField, scalar>(this->phiName_);
172 
173  const fvPatchField<scalar>& psip =
174  patch().lookupPatchField<volScalarField, scalar>(psiName_);
175 
176  scalar gM1ByG = (gamma_ - 1.0)/gamma_;
177 
178  this->refValue() =
179  T0_/(1.0 + 0.5*psip*gM1ByG*(1.0 - pos0(phip))*magSqr(Up));
180  this->valueFraction() = 1.0 - pos0(phip);
181 
182  inletOutletFvPatchScalarField::updateCoeffs();
183 }
184 
185 
187 const
188 {
190  writeEntryIfDifferent<word>(os, "U", "U", UName_);
191  writeEntryIfDifferent<word>(os, "phi", "phi", this->phiName_);
192  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
193  os.writeKeyword("gamma") << gamma_ << token::END_STATEMENT << nl;
194  T0_.writeEntry("T0", os);
195  writeEntry("value", os);
196 }
197 
198 
199 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
200 
201 namespace Foam
202 {
204  (
207  );
208 }
209 
210 // ************************************************************************* //
#define readScalar
Definition: doubleScalar.C:38
Foam::surfaceFields.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
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
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:726
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Macros for easy insertion into run-time selection tables.
This boundary condition provides an outflow condition for total temperature for use with supersonic c...
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:262
inletOutletTotalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
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
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.