Loading...
Searching...
No Matches
createFields.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29 \*---------------------------------------------------------------------------*/
30Info << "Reading thermophysical properties\n" << endl;
31
32IOdictionary transportProperties
33(
34 IOobject
35 (
36 "transportProperties",
37 runTime.constant(),
38 mesh,
39 IOobject::MUST_READ_IF_MODIFIED,
40 IOobject::NO_WRITE
41 )
42);
43
44// Laminar viscocity [m2/s]
45_nu = autoPtr<dimensionedScalar>
46 (
47 new dimensionedScalar
48 (
49 "nu",
50 dimViscosity,
51 transportProperties.lookup("nu")
52 )
53 );
54dimensionedScalar& nu = _nu();
55
56// Turbulent viscocity [m2/s]
57_nut = autoPtr<volScalarField>
58 (
59 new volScalarField
60 (
61 IOobject
62 (
63 "nut",
64 runTime.timeName(),
65 mesh,
66 IOobject::MUST_READ,
67 IOobject::AUTO_WRITE
68 ),
69 mesh
70 )
71 );
72volScalarField& nut = _nut();
73
74_UliftBC = autoPtr<volScalarField>
75 (
76 new volScalarField
77 (
78 IOobject
79 (
80 "UliftBC",
81 runTime.timeName(),
82 mesh,
83 IOobject::MUST_READ,
84 IOobject::AUTO_WRITE
85 ),
86 mesh
87 )
88 );
89volScalarField& UliftBC = _UliftBC();
90
91// Laminar Prandtl number
92_Pr = autoPtr<dimensionedScalar>
93 (
94 new dimensionedScalar
95 (
96 "Pr",
97 dimless,
98 transportProperties.lookup("Pr")
99 )
100 );
101dimensionedScalar& Pr = _Pr();
102
103// Turbulent Prandtl number
104_Prt = autoPtr<dimensionedScalar>
105 (
106 new dimensionedScalar
107 (
108 "Prt",
109 dimless,
110 transportProperties.lookup("Prt")
111 )
112 );
113dimensionedScalar& Prt = _Prt();
114
115// Thermal expansion coefficient [1/K]
116_beta = autoPtr<dimensionedScalar>
117 (
118 new dimensionedScalar
119 (
120 "beta",
121 dimless / dimTemperature,
122 transportProperties.lookup("beta")
123 )
124 );
125dimensionedScalar& beta = _beta();
126
127// Reference temperature [K]
128_TRef = autoPtr<dimensionedScalar>
129 (
130 new dimensionedScalar
131 (
132 "TRef",
133 dimTemperature,
134 transportProperties.lookup("TRef")
135 )
136 );
137dimensionedScalar& TRef = _TRef();
138
139// Gravitational constant
140Info << "\nReading g" << endl;
141_g = autoPtr<dimensionedVector>
142 (
143 new dimensionedVector
144 (
145 "g",
146 dimensionSet(0, 1, -2, 0, 0),
147 vector(0, -9.81, 0)
148 )
149 );
150dimensionedVector& g = _g();
151
152Info << "\nReading hRef" << endl;
153_hRef = autoPtr<dimensionedScalar>
154 (
155 new dimensionedScalar
156 (
157 "hRef",
158 dimLength,
159 0
160 )
161 );
162dimensionedScalar& hRef = _hRef();
163
164Info << "Calculating field ghRef\n" << endl;
165_ghRef = autoPtr<dimensionedScalar>
166 (
167 new dimensionedScalar
168 (
169 "ghRef",
170 dimAcceleration* dimLength,
171 0
172 )
173 );
174dimensionedScalar& ghRef = _ghRef();
175
176Info << "Reading field gh\n" << endl;
177_gh = autoPtr<volScalarField>
178 (
179 new volScalarField
180 (
181 IOobject
182 (
183 "gh",
184 runTime.timeName(),
185 mesh
186 ),
187 (g& mesh.C()) - ghRef
188 )
189 );
190volScalarField& gh = _gh();
191
192Info << "Reading field ghf\n" << endl;
193_ghf = autoPtr<surfaceScalarField>
194 (
195 new surfaceScalarField
196 (
197 IOobject
198 (
199 "ghf",
200 runTime.timeName(),
201 mesh
202 ),
203 (g& mesh.Cf()) - ghRef)
204 );
205surfaceScalarField& ghf = _ghf();
206
207// Reading fields
208Info << "Reading field T\n" << endl;
209_T = autoPtr<volScalarField>
210 (
211 new volScalarField
212 (
213 IOobject
214 (
215 "T",
216 runTime.timeName(),
217 mesh,
218 IOobject::MUST_READ,
219 IOobject::AUTO_WRITE
220 ),
221 mesh
222 )
223 );
224volScalarField& T = _T();
225
226Info << "Reading field p_rgh\n" << endl;
227_p_rgh = autoPtr<volScalarField>
228 (
229 new volScalarField
230 (
231 IOobject
232 (
233 "p_rgh",
234 runTime.timeName(),
235 mesh,
236 IOobject::MUST_READ,
237 IOobject::AUTO_WRITE
238 ),
239 mesh
240 )
241 );
242volScalarField& p_rgh = _p_rgh();
243
244Info << "Reading field U\n" << endl;
245_U = autoPtr<volVectorField>
246 (
247 new volVectorField
248 (
249 IOobject
250 (
251 "U",
252 runTime.timeName(),
253 mesh,
254 IOobject::MUST_READ,
255 IOobject::AUTO_WRITE
256 ),
257 mesh
258 )
259 );
260volVectorField& U = _U();
261
262
263#include "createPhi.H"
264
265_laminarTransport = autoPtr<singlePhaseTransportModel>
266 (
267 new singlePhaseTransportModel( U, phi )
268 );
269singlePhaseTransportModel& laminarTransport = _laminarTransport();
270
271Info << "Creating turbulence model\n" << endl;
272turbulence = autoPtr<incompressible::turbulenceModel>
273 (
274 incompressible::turbulenceModel::New(U, phi, laminarTransport)
275 );
276
277// Kinematic density for buoyancy force
278_rhok = autoPtr<volScalarField>
279 (
280 new volScalarField
281 (
282 IOobject
283 (
284 "rhok",
285 runTime.timeName(),
286 mesh
287 ),
288 1.0 - beta * (T - TRef)
289 )
290 );
291volScalarField& rhok = _rhok();
292
293// kinematic turbulent thermal conductivity m2/s
294Info << "Reading field alphat\n" << endl;
295_alphat = autoPtr<volScalarField>
296 (
297 new volScalarField
298 (
299 IOobject
300 (
301 "alphat",
302 runTime.timeName(),
303 mesh,
304 IOobject::MUST_READ,
305 IOobject::AUTO_WRITE
306 ),
307 mesh
308 )
309 );
310volScalarField& alphat = _alphat();
311
312Info << "Reading field p\n" << endl;
313_p = autoPtr<volScalarField>
314 (
315 new volScalarField
316 (
317 IOobject
318 (
319 "p",
320 runTime.timeName(),
321 mesh,
322 IOobject::NO_READ,
323 IOobject::AUTO_WRITE
324 ),
325 p_rgh + rhok* gh
326 )
327 );
328volScalarField& p = _p();
329
330
333
334pimpleControl& pimple = _pimple();
335
337(
338 p,
339 p_rgh,
340 pimple.dict(),
341 pRefCell,
343);
344
345if (p_rgh.needReference())
346{
347 p += dimensionedScalar
348 (
349 "p",
350 p.dimensions(),
351 pRefValue - getRefCellValue(p, pRefCell)
352 );
353}
354
355mesh.setFluxRequired(p_rgh.name());
356
357_MRF = autoPtr<IOMRFZoneList>
358 (
359 new IOMRFZoneList(mesh)
360 );
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
_T
Definition createT.H:30
surfaceScalarField & phi
volVectorField & U
_nut
volScalarField & gh
volScalarField & p_rgh
pimpleControl & pimple
_beta
_rhok
dimensionedScalar & nu
_TRef
volScalarField & rhok
_UliftBC
_Pr
dimensionedScalar & hRef
dimensionedScalar & ghRef
dimensionedScalar & Prt
volScalarField & nut
volScalarField & UliftBC
dimensionedVector & g
dimensionedScalar & Pr
volScalarField & alphat
dimensionedScalar & TRef
surfaceScalarField & ghf
_Prt
dimensionedScalar & beta
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE))
_nu
volScalarField & T
_laminarTransport
volScalarField & p
scalar pRefValue
turbulence
singlePhaseTransportModel & laminarTransport
_MRF
label pRefCell
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)