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)
44{
45 refValue() = 0.0;
46 refGrad() = 0.0;
47 valueFraction() = 0.0;
48}
49
50
52(
54 const fvPatch& p,
56 const fvPatchFieldMapper& mapper
57)
58:
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() = 0.0;
78 valueFraction() = 0.0;
79
80 if (dict.found("value"))
81 {
83 (
84 scalarField("value", dict, p.size())
85 );
86 }
87 else
88 {
89 evaluate();
90 }
91}
92
93
95(
97)
98:
99 mixedFvPatchScalarField(tppsf),
100 Tinf_(tppsf.Tinf_),
101 alphaWall_(tppsf.alphaWall_)
102{}
103
104
106(
109)
110:
111 mixedFvPatchScalarField(tppsf, iF),
112 Tinf_(tppsf.Tinf_),
113 alphaWall_(tppsf.alphaWall_)
114{}
115
116
117// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
118
120(
121 const fvPatchFieldMapper& m
122)
123{
125 Tinf_.autoMap(m);
126 alphaWall_.autoMap(m);
127}
128
129
131(
132 const fvPatchScalarField& ptf,
133 const labelList& addr
134)
135{
136 mixedFvPatchScalarField::rmap(ptf, addr);
137
139 refCast<const wallHeatTransferFvPatchScalarField>(ptf);
140
141 Tinf_.rmap(tiptf.Tinf_, addr);
142 alphaWall_.rmap(tiptf.alphaWall_, addr);
143}
144
145
147{
148 if (updated())
149 {
150 return;
151 }
152
153 const compressible::turbulenceModel& turbModel =
154 db().lookupObject<compressible::turbulenceModel>
155 (
157 (
159 internalField().group()
160 )
161 );
162
163 const label patchi = patch().index();
164
165 valueFraction() =
166 1.0/
167 (
168 1.0
169 + turbModel.kappaEff(patchi)*patch().deltaCoeffs()/alphaWall_
170 );
171
172 mixedFvPatchScalarField::updateCoeffs();
173}
174
175
177{
179 Tinf_.writeEntry("Tinf", os);
180 alphaWall_.writeEntry("alphaWall", os);
181 writeEntry("value", os);
182}
183
184
185// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
186
187namespace Foam
188{
190 (
193 );
194}
195
196// ************************************************************************* //
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)
Map from self.
Definition: Field.C:403
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:62
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
virtual bool write()
Write the output fields.
A FieldMapper for finite-volume patch fields.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:82
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.
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)
Definition: fvPatchField.H:676
Namespace for OpenFOAM.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dictionary dict