icoFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 Description
28  Transient solver for incompressible, laminar flow of Newtonian fluids.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "pisoControl.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 int main(int argc, char *argv[])
38 {
39  #include "setRootCaseLists.H"
40  #include "createTime.H"
41  #include "createMesh.H"
42 
43  pisoControl piso(mesh);
44 
45  #include "createFields.H"
46  #include "initContinuityErrs.H"
47 
48  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
49 
50  Info<< "\nStarting time loop\n" << endl;
51 
52  while (runTime.loop())
53  {
54  Info<< "Time = " << runTime.timeName() << nl << endl;
55 
56  #include "CourantNo.H"
57 
58  // Momentum predictor
59 
61  (
62  fvm::ddt(U)
63  + fvm::div(phi, U)
64  - fvm::laplacian(nu, U)
65  );
66 
67  if (piso.momentumPredictor())
68  {
69  solve(UEqn == -fvc::grad(p));
70  }
71 
72  // --- PISO loop
73  while (piso.correct())
74  {
75  volScalarField rAU(1.0/UEqn.A());
78  (
79  "phiHbyA",
82  );
83 
84  adjustPhi(phiHbyA, U, p);
85 
86  // Update the pressure BCs to ensure flux consistency
88 
89  // Non-orthogonal pressure corrector loop
90  while (piso.correctNonOrthogonal())
91  {
92  // Pressure corrector
93 
94  fvScalarMatrix pEqn
95  (
97  );
98 
99  pEqn.setReference(pRefCell, pRefValue);
100 
101  pEqn.solve(mesh.solver(p.select(piso.finalInnerIter())));
102 
103  if (piso.finalNonOrthogonalIter())
104  {
105  phi = phiHbyA - pEqn.flux();
106  }
107  }
108 
109  #include "continuityErrs.H"
110 
111  U = HbyA - rAU*fvc::grad(p);
112  U.correctBoundaryConditions();
113  }
114 
115  runTime.write();
116 
117  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
118  << " ClockTime = " << runTime.elapsedClockTime() << " s"
119  << nl << endl;
120  }
121 
122  Info<< "End\n" << endl;
123 
124  return 0;
125 }
126 
127 
128 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
surfaceScalarField & phi
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
void constrainPressure(volScalarField &p, const RhoType &rho, const volVectorField &U, const surfaceScalarField &phiHbyA, const RAUType &rhorAU, const MRFType &MRF)
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
pisoControl piso(mesh)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
tmp< volScalarField > rAU
Definition: initCorrectPhi.H:1
HbyA
Definition: pcEqn.H:74
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:34
phiHbyA
Definition: pcEqn.H:73
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
rhoEqn solve()
tmp< volVectorField > constrainHbyA(const tmp< volVectorField > &tHbyA, const volVectorField &U, const volScalarField &p)
Definition: constrainHbyA.C:33
static const char nl
Definition: Ostream.H:265
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
U
Definition: pEqn.H:72
label pRefCell
Definition: createFields.H:106
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
volScalarField & nu
scalar pRefValue
Definition: createFields.H:107