pimpleFoam.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-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
27 Application
28  pimpleFoam.C
29 
30 Group
31  grpIncompressibleSolvers
32 
33 Description
34  Transient solver for incompressible, turbulent flow of Newtonian fluids
35  on a moving mesh.
36 
37  \heading Solver details
38  The solver uses the PIMPLE (merged PISO-SIMPLE) algorithm to solve the
39  continuity equation:
40 
41  \f[
42  \div \vec{U} = 0
43  \f]
44 
45  and momentum equation:
46 
47  \f[
48  \ddt{\vec{U}} + \div \left( \vec{U} \vec{U} \right) - \div \gvec{R}
49  = - \grad p + \vec{S}_U
50  \f]
51 
52  Where:
53  \vartable
54  \vec{U} | Velocity
55  p | Pressure
56  \vec{R} | Stress tensor
57  \vec{S}_U | Momentum source
58  \endvartable
59 
60  Sub-models include:
61  - turbulence modelling, i.e. laminar, RAS or LES
62  - run-time selectable MRF and finite volume options, e.g. explicit porosity
63 
64  \heading Required fields
65  \plaintable
66  U | Velocity [m/s]
67  p | Kinematic pressure, p/rho [m2/s2]
68  <turbulence fields> | As required by user selection
69  \endplaintable
70 
71 Note
72  The motion frequency of this solver can be influenced by the presence
73  of "updateControl" and "updateInterval" in the dynamicMeshDict.
74 
75 \*---------------------------------------------------------------------------*/
76 
77 #include "fvCFD.H"
78 #include "dynamicFvMesh.H"
81 #include "pimpleControl.H"
82 #include "CorrectPhi.H"
83 #include "fvOptions.H"
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 int main(int argc, char *argv[])
88 {
89  argList::addNote
90  (
91  "Transient solver for incompressible, turbulent flow"
92  " of Newtonian fluids on a moving mesh."
93  );
94 
95  #include "postProcess.H"
96 
97  #include "addCheckCaseOptions.H"
98  #include "setRootCaseLists.H"
99  #include "createTime.H"
100  #include "createDynamicFvMesh.H"
101  #include "initContinuityErrs.H"
102  #include "createDyMControls.H"
103  #include "createFields.H"
104  #include "createUfIfPresent.H"
105  #include "CourantNo.H"
106  #include "setInitialDeltaT.H"
107 
108  turbulence->validate();
109 
110  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
111 
112  Info<< "\nStarting time loop\n" << endl;
113 
114  while (runTime.run())
115  {
116  #include "readDyMControls.H"
117  #include "CourantNo.H"
118  #include "setDeltaT.H"
119 
120  ++runTime;
121 
122  Info<< "Time = " << runTime.timeName() << nl << endl;
123 
124  // --- Pressure-velocity PIMPLE corrector loop
125  while (pimple.loop())
126  {
127  if (pimple.firstIter() || moveMeshOuterCorrectors)
128  {
129  // Do any mesh changes
130  mesh.controlledUpdate();
131 
132  if (mesh.changing())
133  {
134  MRF.update();
135 
136  if (correctPhi)
137  {
138  // Calculate absolute flux
139  // from the mapped surface velocity
140  phi = mesh.Sf() & Uf();
141 
142  #include "correctPhi.H"
143 
144  // Make the flux relative to the mesh motion
146  }
147 
148  if (checkMeshCourantNo)
149  {
150  #include "meshCourantNo.H"
151  }
152  }
153  }
154 
155  #include "UEqn.H"
156 
157  // --- Pressure corrector loop
158  while (pimple.correct())
159  {
160  #include "pEqn.H"
161  }
162 
163  if (pimple.turbCorr())
164  {
165  laminarTransport.correct();
166  turbulence->correct();
167  }
168  }
169 
170  runTime.write();
171 
172  runTime.printExecutionTime(Info);
173  }
174 
175  Info<< "End\n" << endl;
176 
177  return 0;
178 }
179 
180 
181 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
fvOptions.H
turbulence
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\n"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
singlePhaseTransportModel.H
turbulentTransportModel.H
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
correctPhi
correctPhi
Definition: readDyMControls.H:3
pimpleControl.H
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
MRF
IOMRFZoneList & MRF
Definition: setRegionFluidFields.H:22
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
setRootCaseLists.H
addCheckCaseOptions.H
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
readDyMControls.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
postProcess.H
Execute application functionObjects to post-process existing results.
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:385
meshCourantNo.H
Calculates and outputs the mean and maximum Courant Numbers.
CorrectPhi.H
createTime.H
createUfIfPresent.H
Creates and initialises the velocity field Uf if required.
dynamicFvMesh.H
fvCFD.H
checkMeshCourantNo
checkMeshCourantNo
Definition: readDyMControls.H:9
moveMeshOuterCorrectors
moveMeshOuterCorrectors
Definition: readDyMControls.H:15
laminarTransport
singlePhaseTransportModel laminarTransport(U, phi)
createDynamicFvMesh.H