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
sonicFoam
sonicDyMFoam
sonicDyMFoam.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
sonicDyMFoam
26
27
Group
28
grpCompressibleSolvers grpMovingMeshSolvers
29
30
Description
31
Transient solver for trans-sonic/supersonic, turbulent flow of a
32
compressible gas, with optional mesh motion and mesh topology changes.
33
34
\*---------------------------------------------------------------------------*/
35
36
#include "
fvCFD.H
"
37
#include "
dynamicFvMesh.H
"
38
#include "
psiThermo.H
"
39
#include "
turbulentFluidThermoModel.H
"
40
#include "
pimpleControl.H
"
41
#include "
CorrectPhi.H
"
42
#include "
fvOptions.H
"
43
44
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
45
46
int
main(
int
argc,
char
*argv[])
47
{
48
#include "
postProcess.H
"
49
50
#include "
setRootCase.H
"
51
#include "
createTime.H
"
52
#include "
createDynamicFvMesh.H
"
53
#include "createControl.H"
54
#include "createControls.H"
55
#include "createFields.H"
56
#include "createFieldRefs.H"
57
#include "
createFvOptions.H
"
58
#include "
createRhoUf.H
"
59
#include "compressibleCourantNo.H"
60
#include "setInitialDeltaT.H"
61
#include "initContinuityErrs.H"
62
63
turbulence
->validate();
64
65
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67
Info
<<
"\nStarting time loop\n"
<<
endl
;
68
69
while
(runTime.run())
70
{
71
#include "readControls.H"
72
73
{
74
// Store divrhoU from the previous mesh so that it can be mapped
75
// and used in correctPhi to ensure the corrected phi has the
76
// same divergence
77
volScalarField
divrhoU
78
(
79
"divrhoU"
,
80
fvc::div
(
fvc::absolute
(
phi
,
rho
,
U
))
81
);
82
83
#include "compressibleCourantNo.H"
84
#include "setDeltaT.H"
85
86
runTime++;
87
88
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
89
90
// Store momentum to set rhoUf for introduced faces.
91
volVectorField
rhoU(
"rhoU"
,
rho
*
U
);
92
93
// Do any mesh changes
94
mesh
.update();
95
96
if
(
mesh
.changing() &&
correctPhi
)
97
{
98
// Calculate absolute flux from the mapped surface velocity
99
phi
=
mesh
.Sf() &
rhoUf
;
100
101
#include "correctPhi.H"
102
103
// Make the fluxes relative to the mesh-motion
104
fvc::makeRelative
(
phi
,
rho
,
U
);
105
}
106
}
107
108
if
(
mesh
.changing() &&
checkMeshCourantNo
)
109
{
110
#include "
meshCourantNo.H
"
111
}
112
113
#include "rhoEqn.H"
114
Info
<<
"rhoEqn max/min : "
<<
max
(
rho
).value()
115
<<
" "
<<
min
(
rho
).value() <<
endl
;
116
117
// --- Pressure-velocity PIMPLE corrector loop
118
while
(
pimple
.loop())
119
{
120
#include "UEqn.H"
121
#include "EEqn.H"
122
123
// --- Pressure corrector loop
124
while
(
pimple
.correct())
125
{
126
#include "pEqn.H"
127
}
128
129
if
(
pimple
.turbCorr())
130
{
131
turbulence
->correct();
132
}
133
}
134
135
runTime.write();
136
137
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
138
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
139
<<
nl
<<
endl
;
140
}
141
142
Info
<<
"End\n"
<<
endl
;
143
144
return
0;
145
}
146
147
148
// ************************************************************************* //
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
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.
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
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
createRhoUf.H
Creates and initialises the velocity velocity field rhoUf.
Generated by
1.8.11