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 | Copyright (C) 2011-2015 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  icoFoam
26 
27 Group
28  grpIncompressibleSolvers
29 
30 Description
31  Transient solver for incompressible, laminar flow of Newtonian fluids.
32 
33  \heading Solver details
34  The solver uses the PISO algorithm to solve the continuity equation:
35 
36  \f[
37  \div \vec{U} = 0
38  \f]
39 
40  and momentum equation:
41 
42  \f[
43  \ddt{\vec{U}}
44  + \div \left( \vec{U} \vec{U} \right)
45  - \div \left(\nu \grad \vec{U} \right)
46  = - \grad p
47  \f]
48 
49  Where:
50  \vartable
51  \vec{U} | Velocity
52  p | Pressure
53  \endvartable
54 
55  \heading Required fields
56  \plaintable
57  U | Velocity [m/s]
58  p | Kinematic pressure, p/rho [m2/s2]
59  \endplaintable
60 
61 \*---------------------------------------------------------------------------*/
62 
63 #include "fvCFD.H"
64 #include "pisoControl.H"
65 
66 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68 int main(int argc, char *argv[])
69 {
70  #include "setRootCase.H"
71  #include "createTime.H"
72  #include "createMesh.H"
73 
74  pisoControl piso(mesh);
75 
76  #include "createFields.H"
77  #include "initContinuityErrs.H"
78 
79  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81  Info<< "\nStarting time loop\n" << endl;
82 
83  while (runTime.loop())
84  {
85  Info<< "Time = " << runTime.timeName() << nl << endl;
86 
87  #include "CourantNo.H"
88 
89  // Momentum predictor
90 
92  (
93  fvm::ddt(U)
94  + fvm::div(phi, U)
95  - fvm::laplacian(nu, U)
96  );
97 
98  if (piso.momentumPredictor())
99  {
100  solve(UEqn == -fvc::grad(p));
101  }
102 
103  // --- PISO loop
104  while (piso.correct())
105  {
106  volScalarField rAU(1.0/UEqn.A());
107 
108  volVectorField HbyA("HbyA", U);
109  HbyA = rAU*UEqn.H();
111  (
112  "phiHbyA",
113  (fvc::interpolate(HbyA) & mesh.Sf())
115  );
116 
117  adjustPhi(phiHbyA, U, p);
118 
119  // Non-orthogonal pressure corrector loop
120  while (piso.correctNonOrthogonal())
121  {
122  // Pressure corrector
123 
124  fvScalarMatrix pEqn
125  (
127  );
128 
129  pEqn.setReference(pRefCell, pRefValue);
130 
131  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
132 
133  if (piso.finalNonOrthogonalIter())
134  {
135  phi = phiHbyA - pEqn.flux();
136  }
137  }
138 
139  #include "continuityErrs.H"
140 
141  U = HbyA - rAU*fvc::grad(p);
142  U.correctBoundaryConditions();
143  }
144 
145  runTime.write();
146 
147  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
148  << " ClockTime = " << runTime.elapsedClockTime() << " s"
149  << nl << endl;
150  }
151 
152  Info<< "End\n" << endl;
153 
154  return 0;
155 }
156 
157 
158 // ************************************************************************* //
pRefCell
const label pRefCell
Definition: setRegionFluidFields.H:27
pisoControl.H
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
U
U
Definition: pEqn.H:46
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::fvVectorMatrix
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::fvc::ddtCorr
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
rAU
volScalarField rAU("rAU", 1.0/UEqn().A())
setRootCase.H
HbyA
HbyA
Definition: pEqn.H:4
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createMesh.H
createTime.H
fvCFD.H
phiHbyA
phiHbyA
Definition: pEqn.H:21
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
pRefValue
const scalar pRefValue
Definition: setRegionFluidFields.H:28