Loading...
Searching...
No Matches
icoFoam.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-2016 OpenFOAM Foundation
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
27Application
28 icoFoam
29
30Group
31 grpIncompressibleSolvers
32
33Description
34 Transient solver for incompressible, laminar flow of Newtonian fluids.
35
36 \heading Solver details
37 The solver uses the PISO algorithm to solve the continuity equation:
38
39 \f[
40 \div \vec{U} = 0
41 \f]
42
43 and momentum equation:
44
45 \f[
46 \ddt{\vec{U}}
47 + \div \left( \vec{U} \vec{U} \right)
48 - \div \left(\nu \grad \vec{U} \right)
49 = - \grad p
50 \f]
51
52 Where:
53 \vartable
54 \vec{U} | Velocity
55 p | Pressure
56 \endvartable
57
58 \heading Required fields
59 \plaintable
60 U | Velocity [m/s]
61 p | Kinematic pressure, p/rho [m2/s2]
62 \endplaintable
63
64\*---------------------------------------------------------------------------*/
65
66#include "fvCFD.H"
67#include "pisoControl.H"
68
69// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70
71int main(int argc, char *argv[])
72{
73 argList::addNote
74 (
75 "Transient solver for incompressible, laminar flow"
76 " of Newtonian fluids."
77 );
78
79 #include "postProcess.H"
80
81 #include "addCheckCaseOptions.H"
82 #include "setRootCaseLists.H"
83 #include "createTime.H"
84 #include "createMesh.H"
85
86 pisoControl piso(mesh);
87
88 #include "createFields.H"
89 #include "initContinuityErrs.H"
90
91 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
92
93 Info<< "\nStarting time loop\n" << endl;
94
95 while (runTime.loop())
96 {
97 Info<< "Time = " << runTime.timeName() << nl << endl;
98
99 #include "CourantNo.H"
100
101 // Momentum predictor
102
104 (
105 fvm::ddt(U)
106 + fvm::div(phi, U)
107 - fvm::laplacian(nu, U)
108 );
109
110 if (piso.momentumPredictor())
111 {
112 solve(UEqn == -fvc::grad(p));
113 }
114
115 // --- PISO loop
116 while (piso.correct())
117 {
118 volScalarField rAU(1.0/UEqn.A());
121 (
122 "phiHbyA",
123 fvc::flux(HbyA)
124 + fvc::interpolate(rAU)*fvc::ddtCorr(U, phi)
125 );
126
127 adjustPhi(phiHbyA, U, p);
128
129 // Update the pressure BCs to ensure flux consistency
131
132 // Non-orthogonal pressure corrector loop
133 while (piso.correctNonOrthogonal())
134 {
135 // Pressure corrector
136
137 fvScalarMatrix pEqn
138 (
139 fvm::laplacian(rAU, p) == fvc::div(phiHbyA)
140 );
141
142 pEqn.setReference(pRefCell, pRefValue);
143
144 pEqn.solve(p.select(piso.finalInnerIter()));
145
146 if (piso.finalNonOrthogonalIter())
147 {
148 phi = phiHbyA - pEqn.flux();
149 }
150 }
151
152 #include "continuityErrs.H"
153
154 U = HbyA - rAU*fvc::grad(p);
155 U.correctBoundaryConditions();
156 }
157
158 runTime.write();
159
160 runTime.printExecutionTime(Info);
161 }
162
163 Info<< "End\n" << endl;
164
165 return 0;
166}
167
168
169// ************************************************************************* //
Required Classes.
const scalar pRefValue
const label pRefCell
U
Definition pEqn.H:72
volScalarField & p
constrainPressure(p_rgh, rho, U, phiHbyA, rhorAUf, MRF)
fvVectorMatrix & UEqn
Definition UEqn.H:13
phiHbyA
Definition pcEqn.H:73
HbyA
Definition pcEqn.H:74
dynamicFvMesh & mesh
engineTime & runTime
Required Classes.
pisoControl piso(mesh)
adjustPhi(phiHbyA, U, p_rgh)
tmp< volScalarField > rAU
GeometricField< vector, fvPatchField, volMesh > volVectorField
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
fvMatrix< vector > fvVectorMatrix
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Execute application functionObjects to post-process existing results.
volScalarField & nu
CEqn solve()