OpenFOAM
4.1
The OpenFOAM Foundation
Main Page
Related Pages
Modules
Namespaces
Classes
Files
File List
File Members
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Modules
Pages
applications
solvers
compressible
rhoPimpleFoam
rhoPimpleDyMFoam
rhoPimpleDyMFoam.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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6
\\/ M anipulation |
7
-------------------------------------------------------------------------------
8
License
9
This file is part of OpenFOAM.
10
11
OpenFOAM is free software: you can redistribute it and/or modify it
12
under the terms of the GNU General Public License as published by
13
the Free Software Foundation, either version 3 of the License, or
14
(at your option) any later version.
15
16
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19
for more details.
20
21
You should have received a copy of the GNU General Public License
22
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23
24
Application
25
rhoPimpleFoam
26
27
Group
28
grpCompressibleSolvers grpMovingMeshSolvers
29
30
Description
31
Transient solver for turbulent flow of compressible fluids for HVAC and
32
similar applications, with optional mesh motion and mesh topology changes.
33
34
Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
35
pseudo-transient simulations.
36
37
\*---------------------------------------------------------------------------*/
38
39
#include "
fvCFD.H
"
40
#include "
dynamicFvMesh.H
"
41
#include "
psiThermo.H
"
42
#include "
turbulentFluidThermoModel.H
"
43
#include "
bound.H
"
44
#include "
pimpleControl.H
"
45
#include "
CorrectPhi.H
"
46
#include "
fvOptions.H
"
47
#include "
localEulerDdtScheme.H
"
48
#include "
fvcSmooth.H
"
49
50
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51
52
int
main(
int
argc,
char
*argv[])
53
{
54
#include "
postProcess.H
"
55
56
#include "
setRootCase.H
"
57
#include "
createTime.H
"
58
#include "
createDynamicFvMesh.H
"
59
#include "createControl.H"
60
#include "
createRDeltaT.H
"
61
#include "initContinuityErrs.H"
62
#include "createFields.H"
63
#include "createFieldRefs.H"
64
#include "
createFvOptions.H
"
65
#include "
createRhoUf.H
"
66
#include "createControls.H"
67
68
turbulence
->validate();
69
70
if
(!
LTS
)
71
{
72
#include "compressibleCourantNo.H"
73
#include "setInitialDeltaT.H"
74
}
75
76
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
77
78
Info
<<
"\nStarting time loop\n"
<<
endl
;
79
80
while
(runTime.run())
81
{
82
#include "readControls.H"
83
84
{
85
// Store divrhoU from the previous mesh so that it can be mapped
86
// and used in correctPhi to ensure the corrected phi has the
87
// same divergence
88
volScalarField
divrhoU
89
(
90
"divrhoU"
,
91
fvc::div
(
fvc::absolute
(
phi
,
rho
,
U
))
92
);
93
94
if
(
LTS
)
95
{
96
#include "setRDeltaT.H"
97
}
98
else
99
{
100
#include "compressibleCourantNo.H"
101
#include "setDeltaT.H"
102
}
103
104
runTime++;
105
106
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
107
108
// Store momentum to set rhoUf for introduced faces.
109
volVectorField
rhoU(
"rhoU"
,
rho
*
U
);
110
111
// Do any mesh changes
112
mesh
.update();
113
114
if
(
mesh
.changing() &&
correctPhi
)
115
{
116
// Calculate absolute flux from the mapped surface velocity
117
phi
=
mesh
.Sf() &
rhoUf
;
118
119
#include "correctPhi.H"
120
121
// Make the fluxes relative to the mesh-motion
122
fvc::makeRelative
(
phi
,
rho
,
U
);
123
}
124
}
125
126
if
(
mesh
.changing() &&
checkMeshCourantNo
)
127
{
128
#include "
meshCourantNo.H
"
129
}
130
131
#include "rhoEqn.H"
132
Info
<<
"rhoEqn max/min : "
<<
max
(
rho
).value()
133
<<
" "
<<
min
(
rho
).value() <<
endl
;
134
135
// --- Pressure-velocity PIMPLE corrector loop
136
while
(
pimple
.loop())
137
{
138
#include "UEqn.H"
139
#include "EEqn.H"
140
141
// --- Pressure corrector loop
142
while
(
pimple
.correct())
143
{
144
#include "pEqn.H"
145
}
146
147
if
(
pimple
.turbCorr())
148
{
149
turbulence
->correct();
150
}
151
}
152
153
runTime.write();
154
155
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
156
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
157
<<
nl
<<
endl
;
158
}
159
160
Info
<<
"End\n"
<<
endl
;
161
162
return
0;
163
}
164
165
166
// ************************************************************************* //
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition:
createFields.H:23
fvCFD.H
phi
surfaceScalarField & phi
Definition:
setRegionFluidFields.H:8
U
U
Definition:
pEqn.H:83
createFvOptions.H
psiThermo.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dynamicFvMesh.H
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition:
fvcDiv.C:47
pimple
const dictionary & pimple
Definition:
readFluidMultiRegionPIMPLEControls.H:1
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition:
Ostream.H:253
createDynamicFvMesh.H
CorrectPhi.H
createTime.H
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition:
volFieldsFwd.H:55
turbulentFluidThermoModel.H
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition:
volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition:
createDynamicFvMesh.H:18
rho
rho
Definition:
readInitialConditions.H:96
Foam::nl
static const char nl
Definition:
Ostream.H:262
bound.H
Bound the given scalar field if it has gone unbounded.
LTS
bool LTS
Definition:
createRDeltaT.H:1
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fvOptions.H
setRootCase.H
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
createRDeltaT.H
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition:
fvcMeshPhi.C:188
pimpleControl.H
rhoUf
rhoUf
Definition:
pEqn.H:91
localEulerDdtScheme.H
Foam::Info
messageStream Info
checkMeshCourantNo
checkMeshCourantNo
Definition:
readControls.H:5
correctPhi
correctPhi
Definition:
readControls.H:3
postProcess.H
Execute application functionObjects to post-process existing results.
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition:
fvcMeshPhi.C:75
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
createRhoUf.H
Creates and initialises the velocity velocity field rhoUf.
Generated by
1.8.11