FreeFOAM The Cross-Platform CFD Toolkit
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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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 
29 #include <finiteVolume/volFields.H>
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
35 (
36  const fvPatch& p,
38 )
39 :
40  mixedFvPatchScalarField(p, iF),
41  Tinf_(p.size(), 0.0),
42  alphaWall_(p.size(), 0.0)
43 {
44  refValue() = 0.0;
45  refGrad() = 0.0;
46  valueFraction() = 0.0;
47 }
48 
49 
51 (
53  const fvPatch& p,
55  const fvPatchFieldMapper& mapper
56 )
57 :
58  mixedFvPatchScalarField(ptf, p, iF, mapper),
59  Tinf_(ptf.Tinf_, mapper),
60  alphaWall_(ptf.alphaWall_, mapper)
61 {}
62 
63 
65 (
66  const fvPatch& p,
68  const dictionary& dict
69 )
70 :
71  mixedFvPatchScalarField(p, iF),
72  Tinf_("Tinf", dict, p.size()),
73  alphaWall_("alphaWall", dict, p.size())
74 {
75  refValue() = Tinf_;
76  refGrad() = 0.0;
77  valueFraction() = 0.0;
78 
79  if (dict.found("value"))
80  {
82  (
83  scalarField("value", dict, p.size())
84  );
85  }
86  else
87  {
88  evaluate();
89  }
90 }
91 
92 
94 (
96 )
97 :
98  mixedFvPatchScalarField(tppsf),
99  Tinf_(tppsf.Tinf_),
100  alphaWall_(tppsf.alphaWall_)
101 {}
102 
103 
105 (
108 )
109 :
110  mixedFvPatchScalarField(tppsf, iF),
111  Tinf_(tppsf.Tinf_),
112  alphaWall_(tppsf.alphaWall_)
113 {}
114 
115 
116 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
117 
119 (
120  const fvPatchFieldMapper& m
121 )
122 {
123  scalarField::autoMap(m);
124  Tinf_.autoMap(m);
125  alphaWall_.autoMap(m);
126 }
127 
128 
130 (
131  const fvPatchScalarField& ptf,
132  const labelList& addr
133 )
134 {
135  mixedFvPatchScalarField::rmap(ptf, addr);
136 
138  refCast<const wallHeatTransferFvPatchScalarField>(ptf);
139 
140  Tinf_.rmap(tiptf.Tinf_, addr);
141  alphaWall_.rmap(tiptf.alphaWall_, addr);
142 }
143 
144 
146 {
147  if (updated())
148  {
149  return;
150  }
151 
152  const basicThermo& thermo = db().lookupObject<basicThermo>
153  (
154  "thermophysicalProperties"
155  );
156 
157  const label patchi = patch().index();
158 
159  const scalarField& Tw = thermo.T().boundaryField()[patchi];
160  scalarField Cpw = thermo.Cp(Tw, patchi);
161 
162  valueFraction() =
163  1.0/
164  (
165  1.0
166  + Cpw*thermo.alpha().boundaryField()[patchi]
167  *patch().deltaCoeffs()/alphaWall_
168  );
169 
171 }
172 
173 
175 {
177  Tinf_.writeEntry("Tinf", os);
178  alphaWall_.writeEntry("alphaWall", os);
179  writeEntry("value", os);
180 }
181 
182 
183 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
184 
185 namespace Foam
186 {
188 }
189 
190 // ************************ vim: set sw=4 sts=4 et: ************************ //