OpenFOAM
v2512
The open source CFD toolbox
Loading...
Searching...
No Matches
buoyantPimpleFoam.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 | www.openfoam.com
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
Copyright (C) 2011-2017 OpenFOAM Foundation
9
Copyright (C) 2021 OpenCFD Ltd.
10
-------------------------------------------------------------------------------
11
License
12
This file is part of OpenFOAM.
13
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
19
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22
for more details.
23
24
You should have received a copy of the GNU General Public License
25
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27
Application
28
buoyantPimpleFoam
29
30
Group
31
grpHeatTransferSolvers
32
33
Description
34
Transient solver for buoyant, turbulent flow of compressible fluids
35
for ventilation and heat-transfer, with optional mesh motion
36
and mesh topology changes.
37
38
Turbulence is modelled using a run-time selectable compressible RAS or
39
LES model.
40
41
\*---------------------------------------------------------------------------*/
42
43
#include "
fvCFD.H
"
44
#include "
dynamicFvMesh.H
"
45
#include "
rhoThermo.H
"
46
#include "
turbulentFluidThermoModel.H
"
47
#include "
radiationModel.H
"
48
#include "
CorrectPhi.H
"
49
#include "
fvOptions.H
"
50
#include "
pimpleControl.H
"
51
#include "
pressureControl.H
"
52
#include "
localEulerDdtScheme.H
"
53
#include "
fvcSmooth.H
"
54
55
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
56
57
int
main(
int
argc,
char
*argv[])
58
{
59
argList::addNote
60
(
61
"Transient solver for buoyant, turbulent fluid flow"
62
" of compressible fluids, including radiation,"
63
" with optional mesh motion and mesh topology changes."
64
);
65
66
#include "
postProcess.H
"
67
68
#include "
addCheckCaseOptions.H
"
69
#include "
setRootCaseLists.H
"
70
#include "
createTime.H
"
71
#include "
createDynamicFvMesh.H
"
72
#include "
createDyMControls.H
"
73
#include "initContinuityErrs.H"
74
#include "
createFields.H
"
75
#include "
createFieldRefs.H
"
76
#include "
createRhoUfIfPresent.H
"
77
78
turbulence
->validate();
79
80
if
(!
LTS
)
81
{
82
#include "compressibleCourantNo.H"
83
#include "setInitialDeltaT.H"
84
}
85
86
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
87
88
Info
<<
"\nStarting time loop\n"
<<
endl
;
89
90
while
(
runTime
.run())
91
{
92
#include "
readDyMControls.H
"
93
94
// Store divrhoU from the previous mesh
95
// so that it can be mapped and used in correctPhi
96
// to ensure the corrected phi has the same divergence
97
autoPtr<volScalarField> divrhoU;
98
if
(
correctPhi
)
99
{
100
divrhoU.reset
101
(
102
new
volScalarField
103
(
104
"divrhoU"
,
105
fvc::div(fvc::absolute(
phi
,
rho
,
U
))
106
)
107
);
108
}
109
110
if
(
LTS
)
111
{
112
#include "setRDeltaT.H"
113
}
114
else
115
{
116
#include "compressibleCourantNo.H"
117
#include "setDeltaT.H"
118
}
119
120
++
runTime
;
121
122
Info
<<
"Time = "
<<
runTime
.timeName() <<
nl
<<
endl
;
123
124
// --- Pressure-velocity PIMPLE corrector loop
125
while
(
pimple
.loop())
126
{
127
if
(
pimple
.firstIter() ||
moveMeshOuterCorrectors
)
128
{
129
// Store momentum to set rhoUf for introduced faces.
130
autoPtr<volVectorField> rhoU;
131
if
(
rhoUf
.valid())
132
{
133
rhoU.reset(
new
volVectorField
(
"rhoU"
,
rho
*
U
));
134
}
135
136
// Do any mesh changes
137
mesh
.update();
138
139
if
(
mesh
.changing())
140
{
141
gh
= (
g
&
mesh
.C()) - ghRef;
142
ghf
= (
g
&
mesh
.Cf()) - ghRef;
143
144
MRF
.update();
145
146
if
(
correctPhi
)
147
{
148
// Calculate absolute flux
149
// from the mapped surface velocity
150
phi
=
mesh
.Sf() &
rhoUf
();
151
152
#include "correctPhi.H"
153
154
// Make the fluxes relative to the mesh-motion
155
fvc::makeRelative(
phi
,
rho
,
U
);
156
}
157
158
if
(
checkMeshCourantNo
)
159
{
160
#include "
meshCourantNo.H
"
161
}
162
}
163
}
164
165
if
(
pimple
.firstIter() && !
pimple
.SIMPLErho())
166
{
167
#include "rhoEqn.H"
168
}
169
170
#include "
UEqn.H
"
171
#include "
EEqn.H
"
172
173
// --- Pressure corrector loop
174
while
(
pimple
.correct())
175
{
176
#include "
pEqn.H
"
177
}
178
179
if
(
pimple
.turbCorr())
180
{
181
turbulence
->correct();
182
}
183
}
184
185
rho
=
thermo
.rho();
186
187
runTime
.write();
188
189
runTime
.printExecutionTime(Info);
190
}
191
192
Info
<<
"End\n"
<<
endl
;
193
194
return
0;
195
}
196
197
198
// ************************************************************************* //
CorrectPhi.H
addCheckCaseOptions.H
Required Classes.
g
const uniformDimensionedVectorField & g
Definition
createFluidFields.H:26
ghf
const surfaceScalarField & ghf
Definition
setRegionFluidFields.H:16
MRF
IOMRFZoneList & MRF
Definition
setRegionFluidFields.H:20
gh
const volScalarField & gh
Definition
setRegionFluidFields.H:15
pimple
pimpleControl & pimple
Definition
setRegionFluidFields.H:56
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
U
U
Definition
pEqn.H:72
phi
phi
Definition
correctPhiFaceMask.H:34
createDyMControls.H
createDynamicFvMesh.H
mesh
dynamicFvMesh & mesh
Definition
createDynamicFvMesh.H:6
runTime
engineTime & runTime
Definition
createEngineTime.H:13
LTS
bool LTS
Definition
createRDeltaT.H:1
createRhoUfIfPresent.H
Creates and initialises the velocity field rhoUf if required.
rhoUf
autoPtr< surfaceVectorField > rhoUf
Definition
createRhoUfIfPresent.H:27
createTime.H
dynamicFvMesh.H
turbulence
compressible::turbulenceModel & turbulence
Definition
setRegionFluidFields.H:28
fvCFD.H
fvOptions.H
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
EEqn.H
UEqn.H
createFieldRefs.H
pEqn.H
localEulerDdtScheme.H
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition
volFieldsFwd.H:76
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere).
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition
Ostream.H:519
Foam::nl
constexpr char nl
The newline '\n' character (0x0a).
Definition
Ostream.H:50
pimpleControl.H
postProcess.H
Execute application functionObjects to post-process existing results.
pressureControl.H
radiationModel.H
readDyMControls.H
checkMeshCourantNo
checkMeshCourantNo
Definition
readDyMControls.H:9
correctPhi
correctPhi
Definition
readDyMControls.H:3
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition
readDyMControls.H:15
rho
rho
Definition
readInitialConditions.H:88
rhoThermo.H
setRootCaseLists.H
createFields.H
turbulentFluidThermoModel.H
applications
solvers
heatTransfer
buoyantPimpleFoam
buoyantPimpleFoam.C
Generated by
1.16.1