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
27(
28 IOobject
29 (
30 "transportProperties",
31 runTime.constant(),
32 mesh,
33 IOobject::MUST_READ_IF_MODIFIED,
34 IOobject::NO_WRITE
35 )
36);
37
38Info << "Reading field p\n" << endl;
39_p = autoPtr<volScalarField>
40 (
41 new volScalarField
42 (
43 IOobject
44 (
45 "p",
46 runTime.timeName(),
47 mesh,
48 IOobject::MUST_READ,
49 IOobject::AUTO_WRITE
50 ),
51 mesh
52 )
53 );
54volScalarField& p = _p();
55
56volScalarField p0(p);
57
58_p0 = autoPtr<volScalarField>
59 (
60 new volScalarField(p0)
61 );
62
63Info << "Reading field U\n" << endl;
64_U = autoPtr<volVectorField>
65 (
66 new volVectorField
67 (
68 IOobject
69 (
70 "U",
71 runTime.timeName(),
72 mesh,
73 IOobject::MUST_READ,
74 IOobject::AUTO_WRITE
75 ),
76 mesh
77 )
78 );
79volVectorField& U = _U();
80
81volVectorField U0(U);
82
83_U0 = autoPtr<volVectorField>
84 (
85 new volVectorField(U0)
86 );
87
88// Laminar viscocity [m2/s]
89_nu = autoPtr<dimensionedScalar>
90 (
91 new dimensionedScalar
92 (
93 "nu",
94 dimViscosity,
95 transportProperties.lookup("nu")
96 )
97 );
98dimensionedScalar& nu = _nu();
99
100_dt = autoPtr<dimensionedScalar>
101 (
102 new dimensionedScalar
103 (
104 "dt",
105 dimensionSet(0, 0, 1, 0, 0),
106 1.0
107 )
108 );
109
110_phi = autoPtr<surfaceScalarField>
111 (
112 new surfaceScalarField
113 (
114 IOobject
115 (
116 "phi",
117 runTime.timeName(),
118 mesh,
119 IOobject::READ_IF_PRESENT,
120 IOobject::AUTO_WRITE
121 ),
122 mesh
123 )
124 );
125surfaceScalarField& phi = _phi();
126
127surfaceScalarField phi0(phi);
128
129_phi0 = autoPtr<surfaceScalarField>
130 (
131 new surfaceScalarField(phi0)
132 );
133
136setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue);
137
138_laminarTransport = autoPtr<singlePhaseTransportModel>
139 (
140 new singlePhaseTransportModel( U, phi )
141 );
142singlePhaseTransportModel& laminarTransport = _laminarTransport();
143
144turbulence = autoPtr<incompressible::turbulenceModel>
145 (
146 incompressible::turbulenceModel::New(U, phi, laminarTransport)
147 );
148
149_MRF = autoPtr<IOMRFZoneList>
150 (
151 new IOMRFZoneList(mesh)
152 );
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
surfaceScalarField & phi
surfaceScalarField phi0(phi)
_U0
volVectorField & U
volVectorField U0(U)
_phi0
dimensionedScalar & nu
_phi
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE))
_dt
_nu
_laminarTransport
volScalarField & p
scalar pRefValue
_p0
turbulence
singlePhaseTransportModel & laminarTransport
_MRF
label pRefCell
volScalarField p0(p)
setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue)