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 License
13 This file is part of ITHACA-FV
14 ITHACA-FV is free software: you can redistribute it and/or modify
15 it under the terms of the GNU Lesser General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18 ITHACA-FV is distributed in the hope that it will be useful,
19 but WITHOUT ANY WARRANTY; without even the implied warranty of
20 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21 GNU Lesser General Public License for more details.
22 You should have received a copy of the GNU Lesser General Public License
23 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
24\*---------------------------------------------------------------------------*/
25
26Info << "Reading field p\n" << endl;
27Info << runTime.timeName() << endl;
28_p = autoPtr<volScalarField>
29 (
30 new volScalarField
31 (
32 IOobject
33 (
34 "p",
35 runTime.timeName(),
36 mesh,
37 IOobject::MUST_READ,
38 IOobject::AUTO_WRITE
39 ),
40 mesh
41 )
42 );
43volScalarField& p = _p();
44
45volScalarField p0(p);
46
47_p0 = autoPtr<volScalarField>
48 (
49 new volScalarField(p0)
50 );
51
52Info << "Reading field U\n" << endl;
53_U = autoPtr<volVectorField>
54 (
55 new volVectorField
56 (
57 IOobject
58 (
59 "U",
60 runTime.timeName(),
61 mesh,
62 IOobject::MUST_READ,
63 IOobject::AUTO_WRITE
64 ),
65 mesh
66 )
67 );
68volVectorField& U = _U();
69
70volVectorField U0(U);
71
72_U0 = autoPtr<volVectorField>
73 (
74 new volVectorField(U0)
75 );
76
77#include "createPhi.H"
78
79
83
84_laminarTransport = autoPtr<singlePhaseTransportModel>
85 (
86 new singlePhaseTransportModel( U, phi )
87 );
88singlePhaseTransportModel& laminarTransport = _laminarTransport();
89
90turbulence = autoPtr<incompressible::turbulenceModel>
91 (
92 incompressible::turbulenceModel::New(U, phi, laminarTransport)
93 );
94
95_MRF = autoPtr<IOMRFZoneList>
96 (
97 new IOMRFZoneList(mesh)
98 );
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
simpleControl simple(mesh)
surfaceScalarField & phi
_U0
volVectorField & U
volVectorField U0(U)
_laminarTransport
volScalarField & p
scalar pRefValue
_p0
turbulence
singlePhaseTransportModel & laminarTransport
_MRF
label pRefCell
volScalarField p0(p)
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)