totalTemperatureFvPatchScalarField.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 #include "surfaceFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38  const dictionary& dict
39 )
40 :
41  fixedValueFvPatchScalarField(p, iF, dict, false),
42  UName_(dict.lookupOrDefault<word>("U", "U")),
43  phiName_(dict.lookupOrDefault<word>("phi", "phi")),
44  psiName_(dict.lookupOrDefault<word>("psi", "psi")),
45  gamma_(dict.lookup<scalar>("gamma", dimless)),
46  T0_("T0", dimTemperature, dict, p.size())
47 {
48  if (dict.found("value"))
49  {
51  (
52  scalarField("value", iF.dimensions(), dict, p.size())
53  );
54  }
55  else
56  {
58  }
59 }
60 
61 
63 (
65  const fvPatch& p,
67  const fieldMapper& mapper
68 )
69 :
70  fixedValueFvPatchScalarField(ptf, p, iF, mapper),
71  UName_(ptf.UName_),
72  phiName_(ptf.phiName_),
73  psiName_(ptf.psiName_),
74  gamma_(ptf.gamma_),
75  T0_(mapper(ptf.T0_))
76 {}
77 
78 
80 (
83 )
84 :
85  fixedValueFvPatchScalarField(tppsf, iF),
86  UName_(tppsf.UName_),
87  phiName_(tppsf.phiName_),
88  psiName_(tppsf.psiName_),
89  gamma_(tppsf.gamma_),
90  T0_(tppsf.T0_)
91 {}
92 
93 
94 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
95 
97 (
98  const fvPatchScalarField& ptf,
99  const fieldMapper& mapper
100 )
101 {
102  fixedValueFvPatchScalarField::map(ptf, mapper);
103 
105  refCast<const totalTemperatureFvPatchScalarField>(ptf);
106 
107  mapper(T0_, tiptf.T0_);
108 }
109 
110 
112 (
113  const fvPatchScalarField& ptf
114 )
115 {
116  fixedValueFvPatchScalarField::reset(ptf);
117 
119  refCast<const totalTemperatureFvPatchScalarField>(ptf);
120 
121  T0_.reset(tiptf.T0_);
122 }
123 
124 
126 {
127  if (updated())
128  {
129  return;
130  }
131 
132  const fvPatchVectorField& Up =
133  patch().lookupPatchField<volVectorField, vector>(UName_);
134 
135  const fvsPatchField<scalar>& phip =
136  patch().lookupPatchField<surfaceScalarField, scalar>(phiName_);
137 
138  const fvPatchField<scalar>& psip =
139  patch().lookupPatchField<volScalarField, scalar>(psiName_);
140 
141  scalar gM1ByG = (gamma_ - 1.0)/gamma_;
142 
143  operator==
144  (
145  T0_/(1.0 + 0.5*psip*gM1ByG*neg(phip)*magSqr(Up))
146  );
147 
148  fixedValueFvPatchScalarField::updateCoeffs();
149 }
150 
151 
153 {
155  writeEntryIfDifferent<word>(os, "U", "U", UName_);
156  writeEntryIfDifferent<word>(os, "phi", "phi", phiName_);
157  writeEntryIfDifferent<word>(os, "psi", "psi", psiName_);
158  writeEntry(os, "gamma", gamma_);
159  writeEntry(os, "T0", T0_);
160  writeEntry(os, "value", *this);
161 }
162 
163 
164 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
165 
166 namespace Foam
167 {
169  (
172  );
173 }
174 
175 // ************************************************************************* //
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
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:83
This boundary condition provides a total temperature condition.
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 fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
totalTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &, const dictionary &)
Construct from patch, internal field and dictionary.
A class for handling words, derived from string.
Definition: word.H:62
Namespace for OpenFOAM.
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimTemperature
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:64
dimensionedScalar neg(const dimensionedScalar &ds)
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dictionary dict
volScalarField & p
Foam::surfaceFields.