Loading...
Searching...
No Matches
wallHeatTransferFvPatchScalarField.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-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30#include "fvPatchFieldMapper.H"
32
33// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34
36(
37 const fvPatch& p,
39)
40:
41 mixedFvPatchScalarField(p, iF),
42 Tinf_(p.size(), Zero),
43 alphaWall_(p.size(), Zero)
45 refValue() = Zero;
46 refGrad() = Zero;
47 valueFraction() = 0.0;
48}
49
50
52(
54 const fvPatch& p,
56 const fvPatchFieldMapper& mapper
57)
59 mixedFvPatchScalarField(ptf, p, iF, mapper),
60 Tinf_(ptf.Tinf_, mapper),
61 alphaWall_(ptf.alphaWall_, mapper)
62{}
63
64
66(
67 const fvPatch& p,
69 const dictionary& dict
70)
71:
72 mixedFvPatchScalarField(p, iF),
73 Tinf_("Tinf", dict, p.size()),
74 alphaWall_("alphaWall", dict, p.size())
75{
76 refValue() = Tinf_;
77 refGrad() = Zero;
78 valueFraction() = 0.0;
79
80 if (!this->readValueEntry(dict))
81 {
82 evaluate();
83 }
84}
85
86
88(
89 const wallHeatTransferFvPatchScalarField& tppsf
90)
92 mixedFvPatchScalarField(tppsf),
93 Tinf_(tppsf.Tinf_),
94 alphaWall_(tppsf.alphaWall_)
95{}
96
97
99(
102)
103:
104 mixedFvPatchScalarField(tppsf, iF),
105 Tinf_(tppsf.Tinf_),
106 alphaWall_(tppsf.alphaWall_)
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
113(
114 const fvPatchFieldMapper& m
115)
118 Tinf_.autoMap(m);
119 alphaWall_.autoMap(m);
120}
121
122
124(
125 const fvPatchScalarField& ptf,
126 const labelList& addr
127)
128{
129 mixedFvPatchScalarField::rmap(ptf, addr);
130
133
134 Tinf_.rmap(tiptf.Tinf_, addr);
135 alphaWall_.rmap(tiptf.alphaWall_, addr);
136}
137
138
140{
141 if (updated())
142 {
143 return;
144 }
145
146 const compressible::turbulenceModel& turbModel =
147 db().lookupObject<compressible::turbulenceModel>
148 (
150 (
152 internalField().group()
153 )
154 );
155
156 const label patchi = patch().index();
157
158 valueFraction() =
159 1.0/
160 (
161 1.0
162 + turbModel.kappaEff(patchi)*patch().deltaCoeffs()/alphaWall_
163 );
164
165 mixedFvPatchScalarField::updateCoeffs();
166}
167
168
170{
172 Tinf_.writeEntry("Tinf", os);
173 alphaWall_.writeEntry("alphaWall", os);
176
177
178// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
179
180namespace Foam
181{
183 (
185 wallHeatTransferFvPatchScalarField
186 );
187}
188
189// ************************************************************************* //
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...
void autoMap(const FieldMapper &map, const bool applyFlip=true)
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A FieldMapper for finite-volume patch fields.
virtual void write(Ostream &) const
Write.
void writeValueEntry(Ostream &os) const
Write *this field as a "value" entry.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
static const word propertiesName
Default name of the turbulence properties dictionary.
This boundary condition provides an enthalpy condition for wall heat transfer.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
wallHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
volScalarField & p
OBJstream os(runTime.globalPath()/outputName)
#define makePatchTypeField(PatchTypeField, typePatchTypeField)
Define a concrete fvPatchField type and add to run-time tables Example, (fvPatchScalarField,...
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
List< label > labelList
A List of labels.
Definition List.H:62
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
fvPatchField< scalar > fvPatchScalarField
dictionary dict