QGDsolver
The open source CFD toolbox
Main Page
Modules
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Friends
Macros
Groups
QHDFoam
createFields.H
Go to the documentation of this file.
1
/*---------------------------------------------------------------------------*\
2
========= |
3
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4
\\ / O peration |
5
\\ / A nd | www.openfoam.com
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
Copyright (C) 2011-2016 OpenFOAM Foundation
9
Copyright (C) 2019 OpenCFD Ltd.
10
Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11
-------------------------------------------------------------------------------
12
License
13
This file is part of QGDsolver library, based on OpenFOAM+.
14
OpenFOAM is free software: you can redistribute it and/or modify it
15
under the terms of the GNU 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
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21
for more details.
22
You should have received a copy of the GNU General Public License
23
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24
Global
25
createFields
26
Description
27
Creating fields for calculations
28
SourceFile
29
QHDFoam.C
30
/*---------------------------------------------------------------------------*\
31
Info<< "Reading thermophysical properties\n" << endl;
32
33
autoPtr<rhoQGDThermo> pThermo
34
(
35
rhoQGDThermo::New(mesh)
36
);
37
rhoQGDThermo& thermo = pThermo();
38
thermo.correct();
39
40
volScalarField& e = thermo.he();
41
42
volScalarField& p = thermo.p();
43
const volScalarField& T = thermo.T();
44
const surfaceScalarField& hQGDf = thermo.hQGDf();
45
const surfaceScalarField& tauQGDf = thermo.tauQGDf();
46
47
Info << "Thermo corrected" << endl;
48
49
autoPtr<compressible::turbulenceModel> turbulence;
50
51
volVectorField U
52
(
53
IOobject
54
(
55
"U",
56
runTime.timeName(),
57
mesh,
58
IOobject::MUST_READ,
59
IOobject::AUTO_WRITE
60
),
61
mesh
62
);
63
64
volScalarField rho
65
(
66
IOobject
67
(
68
"rho",
69
runTime.timeName(),
70
mesh,
71
IOobject::NO_READ,
72
IOobject::AUTO_WRITE
73
),
74
thermo.rho()
75
);
76
77
volVectorField W
78
(
79
IOobject
80
(
81
"W",
82
runTime.timeName(),
83
mesh,
84
IOobject::NO_READ,
85
IOobject::NO_WRITE
86
),
87
U
88
);
89
90
91
92
93
Info<< "\nReading gravitationalProperties" << endl;
94
95
IOdictionary gravitationalProperties
96
(
97
IOobject
98
(
99
"gravitationalProperties",
100
runTime.constant(),
101
mesh,
102
IOobject::MUST_READ_IF_MODIFIED,
103
IOobject::NO_WRITE
104
)
105
);
106
107
const dimensionedVector g(gravitationalProperties.lookup("g"));
108
109
dimensionedScalar beta
110
(
111
"beta",
112
dimless/dimTemperature,
113
thermo.subDict("mixture").subDict("transport")
114
);
115
116
117
surfaceScalarField phiu
118
(
119
"phiu",
120
mesh.Sf() & linearInterpolate(U)
121
);
122
123
surfaceScalarField phiwo
124
(
125
"phiwo",
126
mesh.Sf() & linearInterpolate(W)
127
);
128
129
surfaceScalarField phi
130
(
131
"phi",
132
mesh.Sf() & (linearInterpolate(U) - linearInterpolate(W))
133
);
134
135
volVectorField BdFrc
136
(
137
"BdFrc",
138
T*g*beta
139
);
140
141
Switch implicitDiffusion (thermo.implicitDiffusion());
142
143
surfaceScalarField phiRhof
144
(
145
"phiRhof",
146
linearInterpolate(rho)*phi
147
);
148
149
Info<< "Creating turbulence model\n" << endl;
150
turbulence.reset
151
(
152
compressible::turbulenceModel::New
153
(
154
rho,
155
U,
156
phiRhof,
157
thermo
158
).ptr()
159
);
160
161
label pRefCell = 0;
162
scalar pRefValue = 0.0;
163
setRefCell(p, mesh.solutionDict().subDict("PISO"), pRefCell, pRefValue);
164
165
#include "createIcoZeroSources.H"
166
167
//
168
//END-OF-FILE
169
//
170
Generated by
1.8.5