Loading...
Searching...
No Matches
turbulentHeatFluxTemperatureFvPatchScalarField.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-2015 OpenFOAM Foundation
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8License
9 This file is part of OpenFOAM.
10 OpenFOAM is free software: you can redistribute it and/or modify it
11 under the terms of the GNU General Public License as published by
12 the Free Software Foundation, either version 3 of the License, or
13 (at your option) any later version.
14 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
15 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
16 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
17 for more details.
18 You should have received a copy of the GNU General Public License
19 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
20\*---------------------------------------------------------------------------*/
21
23#include "addToRunTimeSelectionTable.H"
24#include "fvPatchFieldMapper.H"
25#include "turbulentTransportModel.H"
26#include "volFields.H"
27#include "IncompressibleTurbulenceModel.H"
28//#include "incompressible/IncomperssibleTurbulenceModel/IncompressibleTurbulenceModel.H"
29
30// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31
32namespace Foam
33{
34// declare specialization within 'Foam' namespace
35template<>
36const char* NamedEnum
37<
40 2
41 >::names[] =
42{
43 "power",
44 "flux"
45};
46}
47
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51
52namespace Foam
53{
54
55namespace incompressible
56{
57
58// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
59
60const NamedEnum
61<
63 2
64 > turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceTypeNames_;
65
66
67// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68
71(
72 const fvPatch& p,
73 const DimensionedField<scalar, volMesh>& iF
74)
75 :
76 fixedGradientFvPatchScalarField(p, iF),
77 heatSource_(hsPower),
78 q_(p.size(), 0.0),
79 alphaEffName_("undefinedAlphaEff")
80{}
81
82
85(
87 const fvPatch& p,
88 const DimensionedField<scalar, volMesh>& iF,
89 const fvPatchFieldMapper& mapper
90)
91 :
92 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
93 heatSource_(ptf.heatSource_),
94 q_(ptf.q_, mapper),
95 alphaEffName_(ptf.alphaEffName_)
96{}
97
98
101(
102 const fvPatch& p,
103 const DimensionedField<scalar, volMesh>& iF,
104 const dictionary& dict
105)
106 :
107 fixedGradientFvPatchScalarField(p, iF),
108 heatSource_(heatSourceTypeNames_.read(dict.lookup("heatSource"))),
109 q_("q", dict, p.size()),
110 alphaEffName_(dict.lookup("alphaEff"))
111{
112 if (dict.found("value") && dict.found("gradient"))
113 {
114 fvPatchField<scalar>::operator=(Field<scalar>("value", dict, p.size()));
115 gradient() = Field<scalar>("gradient", dict, p.size());
116 }
117 else
118 {
119 fvPatchField<scalar>::operator=(patchInternalField());
120 gradient() = 0.0;
121 }
122}
123
124
127(
129)
130 :
131 fixedGradientFvPatchScalarField(thftpsf),
132 heatSource_(thftpsf.heatSource_),
133 q_(thftpsf.q_),
134 alphaEffName_(thftpsf.alphaEffName_)
135{}
136
137
140(
142 const DimensionedField<scalar, volMesh>& iF
143)
144 :
145 fixedGradientFvPatchScalarField(thftpsf, iF),
146 heatSource_(thftpsf.heatSource_),
147 q_(thftpsf.q_),
148 alphaEffName_(thftpsf.alphaEffName_)
149{}
150
151
152// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
153
155(
156 const fvPatchFieldMapper& m
157)
158{
159 scalarField::autoMap(m);
160 q_.autoMap(m);
161}
162
163
165(
166 const fvPatchScalarField& ptf,
167 const labelList& addr
168)
169{
170 fixedGradientFvPatchScalarField::rmap(ptf, addr);
172 refCast<const turbulentHeatFluxTemperatureFvPatchScalarField>
173 (
174 ptf
175 );
176 q_.rmap(thftptf.q_, addr);
177}
178
179
181{
182 if (updated())
183 {
184 return;
185 }
186
187 const scalarField& alphaEffp =
188 patch().lookupPatchField<volScalarField, scalar>(alphaEffName_);
189 // retrieve (constant) specific heat capacity from transport dictionary
190 const IOdictionary& transportProperties =
191 db().lookupObject<IOdictionary>("transportProperties");
192 const scalar rhoCp0(readScalar(transportProperties.lookup("rhoCp0")));
193
194 switch (heatSource_)
195 {
196 case hsPower:
197 {
198 const scalar Ap = gSum(patch().magSf());
199 gradient() = q_ / (Ap * rhoCp0 * alphaEffp);
200 break;
201 }
202
203 case hsFlux:
204 {
205 gradient() = q_ / (rhoCp0 * alphaEffp);
206 break;
207 }
208
209 default:
210 {
211 FatalErrorIn
212 (
213 "turbulentHeatFluxTemperatureFvPatchScalarField"
214 "("
215 "const fvPatch&, "
216 "const DimensionedField<scalar, volMesh>&, "
217 "const dictionary&"
218 ")"
219 ) << "Unknown heat source type. Valid types are: "
220 << heatSourceTypeNames_ << nl << exit(FatalError);
221 }
222 }
223
224 fixedGradientFvPatchScalarField::updateCoeffs();
225}
226
227
229{
230 fixedGradientFvPatchScalarField::write(os);
231 os.writeKeyword("heatSource") << heatSourceTypeNames_[heatSource_]
232 << token::END_STATEMENT << nl;
233 q_.writeEntry("q", os);
234 os.writeKeyword("alphaEff") << alphaEffName_ << token::END_STATEMENT << nl;
235 writeEntry("value", os);
236}
237
238
239// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
240
242(
243 fvPatchScalarField,
245);
246
247
248// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249
250} // End namespace incompressible
251} // End namespace Foam
252
253// ************************************************************************* //
turbulence read()
turbulentHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
makePatchTypeField(fvPatchScalarField, turbulentHeatFluxTemperatureFvPatchScalarField)
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE))
volScalarField & p