wallHeatTransferFvPatchScalarField.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-2019 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"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
34 (
35  const fvPatch& p,
37 )
38 :
39  mixedFvPatchScalarField(p, iF),
40  Tinf_(p.size(), 0.0),
41  alphaWall_(p.size(), 0.0)
42 {
43  refValue() = 0.0;
44  refGrad() = 0.0;
45  valueFraction() = 0.0;
46 }
47 
48 
50 (
52  const fvPatch& p,
54  const fvPatchFieldMapper& mapper
55 )
56 :
57  mixedFvPatchScalarField(ptf, p, iF, mapper),
58  Tinf_(mapper(ptf.Tinf_)),
59  alphaWall_(mapper(ptf.alphaWall_))
60 {}
61 
62 
64 (
65  const fvPatch& p,
67  const dictionary& dict
68 )
69 :
70  mixedFvPatchScalarField(p, iF),
71  Tinf_("Tinf", dict, p.size()),
72  alphaWall_("alphaWall", dict, p.size())
73 {
74  refValue() = Tinf_;
75  refGrad() = 0.0;
76  valueFraction() = 0.0;
77 
78  if (dict.found("value"))
79  {
81  (
82  scalarField("value", dict, p.size())
83  );
84  }
85  else
86  {
87  evaluate();
88  }
89 }
90 
91 
93 (
95 )
96 :
97  mixedFvPatchScalarField(tppsf),
98  Tinf_(tppsf.Tinf_),
99  alphaWall_(tppsf.alphaWall_)
100 {}
101 
102 
104 (
107 )
108 :
109  mixedFvPatchScalarField(tppsf, iF),
110  Tinf_(tppsf.Tinf_),
111  alphaWall_(tppsf.alphaWall_)
112 {}
113 
114 
115 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
116 
118 (
119  const fvPatchFieldMapper& m
120 )
121 {
122  m(*this, *this);
123  m(Tinf_, Tinf_);
124  m(alphaWall_, alphaWall_);
125 }
126 
127 
129 (
130  const fvPatchScalarField& ptf,
131  const labelList& addr
132 )
133 {
134  mixedFvPatchScalarField::rmap(ptf, addr);
135 
137  refCast<const wallHeatTransferFvPatchScalarField>(ptf);
138 
139  Tinf_.rmap(tiptf.Tinf_, addr);
140  alphaWall_.rmap(tiptf.alphaWall_, addr);
141 }
142 
143 
145 {
146  if (updated())
147  {
148  return;
149  }
150 
151  const compressible::turbulenceModel& turbModel =
152  db().lookupObject<compressible::turbulenceModel>
153  (
155  (
157  internalField().group()
158  )
159  );
160 
161  const label patchi = patch().index();
162 
163  valueFraction() =
164  1.0/
165  (
166  1.0
167  + turbModel.kappaEff(patchi)*patch().deltaCoeffs()/alphaWall_
168  );
169 
170  mixedFvPatchScalarField::updateCoeffs();
171 }
172 
173 
175 {
177  writeEntry(os, "Tinf", Tinf_);
178  writeEntry(os, "alphaWall", alphaWall_);
179  writeEntry(os, "value", *this);
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 namespace Foam
186 {
188  (
191  );
192 }
193 
194 // ************************************************************************* //
const char *const group
Group name for atomic constants.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual tmp< volScalarField > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
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
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:280
Macros for easy insertion into run-time selection tables.
static const word propertiesName
Default name of the turbulence properties dictionary.
static word groupName(Name name, const word &group)
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:155
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
volScalarField scalarField(fieldObject, mesh)
This boundary condition provides an enthalpy condition for wall heat transfer.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
label patchi
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
wallHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
Namespace for OpenFOAM.