Home
Downloads
Documentation
Installation
User Guide
man-pages
API Documentation
README
Release Notes
Changes
License
Support
SourceForge Project
Main Page
Related Pages
Namespaces
Classes
Files
File List
File Members
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Friends
Macros
Pages
applications
solvers
incompressible
pisoFoam
pisoFoam.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) 1991-2010 OpenCFD Ltd.
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
pisoFoam
26
27
Description
28
Transient solver for incompressible flow.
29
30
Turbulence modelling is generic, i.e. laminar, RAS or LES may be selected.
31
32
Usage
33
- pisoFoam [OPTION]
34
35
@param -case <dir> \n
36
Specify the case directory
37
38
@param -parallel \n
39
Run the case in parallel
40
41
@param -help \n
42
Display short usage message
43
44
@param -doc \n
45
Display Doxygen documentation page
46
47
@param -srcDoc \n
48
Display source code
49
50
\*---------------------------------------------------------------------------*/
51
52
#include <
finiteVolume/fvCFD.H
>
53
#include <
incompressibleTransportModels/singlePhaseTransportModel.H
>
54
#include <
incompressibleTurbulenceModel/turbulenceModel.H
>
55
56
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
57
58
int
main
(
int
argc,
char
*argv[])
59
{
60
#include <
OpenFOAM/setRootCase.H
>
61
62
#include <
OpenFOAM/createTime.H
>
63
#include <
OpenFOAM/createMesh.H
>
64
#include "
createFields.H
"
65
#include <
finiteVolume/initContinuityErrs.H
>
66
67
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
68
69
Info
<<
"\nStarting time loop\n"
<<
endl
;
70
71
while
(runTime.loop())
72
{
73
Info
<<
"Time = "
<< runTime.timeName() <<
nl
<<
endl
;
74
75
#include <
finiteVolume/readPISOControls.H
>
76
#include <
finiteVolume/CourantNo.H
>
77
78
// Pressure-velocity PISO corrector
79
{
80
// Momentum predictor
81
82
fvVectorMatrix
UEqn
83
(
84
fvm::ddt
(
U
)
85
+
fvm::div
(
phi
,
U
)
86
+
turbulence
->divDevReff(
U
)
87
);
88
89
UEqn
.
relax
();
90
91
if
(
momentumPredictor
)
92
{
93
solve
(
UEqn
== -
fvc::grad
(
p
));
94
}
95
96
// --- PISO loop
97
98
for
(
int
corr=0; corr<
nCorr
; corr++)
99
{
100
volScalarField
rUA
= 1.0/
UEqn
.
A
();
101
102
U
= rUA*
UEqn
.
H
();
103
phi
= (
fvc::interpolate
(
U
) &
mesh
.
Sf
())
104
+
fvc::ddtPhiCorr
(rUA,
U
,
phi
);
105
106
adjustPhi
(
phi
,
U
,
p
);
107
108
// Non-orthogonal pressure corrector loop
109
for
(
int
nonOrth=0; nonOrth<=
nNonOrthCorr
; nonOrth++)
110
{
111
// Pressure corrector
112
113
fvScalarMatrix
pEqn
114
(
115
fvm::laplacian
(rUA,
p
) ==
fvc::div
(
phi
)
116
);
117
118
pEqn.setReference(
pRefCell
,
pRefValue
);
119
120
if
121
(
122
corr == nCorr-1
123
&& nonOrth ==
nNonOrthCorr
124
)
125
{
126
pEqn.solve(
mesh
.
solver
(
"pFinal"
));
127
}
128
else
129
{
130
pEqn.solve();
131
}
132
133
if
(nonOrth ==
nNonOrthCorr
)
134
{
135
phi
-= pEqn.flux();
136
}
137
}
138
139
#include <
finiteVolume/continuityErrs.H
>
140
141
U
-= rUA*
fvc::grad
(
p
);
142
U
.correctBoundaryConditions();
143
}
144
}
145
146
turbulence
->correct();
147
148
runTime.write();
149
150
Info
<<
"ExecutionTime = "
<< runTime.elapsedCpuTime() <<
" s"
151
<<
" ClockTime = "
<< runTime.elapsedClockTime() <<
" s"
152
<<
nl
<<
endl
;
153
}
154
155
Info
<<
"End\n"
<<
endl
;
156
157
return
0;
158
}
159
160
161
// ************************ vim: set sw=4 sts=4 et: ************************ //