rhoPimpleFoam.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-2016 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  rhoPimpleFoam
26 
27 Description
28  Transient solver for turbulent flow of compressible fluids for HVAC and
29  similar applications.
30 
31  Uses the flexible PIMPLE (PISO-SIMPLE) solution for time-resolved and
32  pseudo-transient simulations.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "psiThermo.H"
39 #include "bound.H"
40 #include "pimpleControl.H"
41 #include "fvOptions.H"
42 #include "localEulerDdtScheme.H"
43 #include "fvcSmooth.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "postProcess.H"
50 
51  #include "setRootCase.H"
52  #include "createTime.H"
53  #include "createMesh.H"
54  #include "createControl.H"
55  #include "createTimeControls.H"
56  #include "createRDeltaT.H"
57  #include "initContinuityErrs.H"
58  #include "createFields.H"
59  #include "createFieldRefs.H"
60  #include "createFvOptions.H"
61 
62  turbulence->validate();
63 
64  if (!LTS)
65  {
66  #include "compressibleCourantNo.H"
67  #include "setInitialDeltaT.H"
68  }
69 
70  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
71 
72  Info<< "\nStarting time loop\n" << endl;
73 
74  while (runTime.run())
75  {
76  #include "readTimeControls.H"
77 
78  if (LTS)
79  {
80  #include "setRDeltaT.H"
81  }
82  else
83  {
84  #include "compressibleCourantNo.H"
85  #include "setDeltaT.H"
86  }
87 
88  runTime++;
89 
90  Info<< "Time = " << runTime.timeName() << nl << endl;
91 
92  if (pimple.nCorrPIMPLE() <= 1)
93  {
94  #include "rhoEqn.H"
95  }
96 
97  // --- Pressure-velocity PIMPLE corrector loop
98  while (pimple.loop())
99  {
100  #include "UEqn.H"
101  #include "EEqn.H"
102 
103  // --- Pressure corrector loop
104  while (pimple.correct())
105  {
106  if (pimple.consistent())
107  {
108  #include "pcEqn.H"
109  }
110  else
111  {
112  #include "pEqn.H"
113  }
114  }
115 
116  if (pimple.turbCorr())
117  {
118  turbulence->correct();
119  }
120  }
121 
122  runTime.write();
123 
124  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
125  << " ClockTime = " << runTime.elapsedClockTime() << " s"
126  << nl << endl;
127  }
128 
129  Info<< "End\n" << endl;
130 
131  return 0;
132 }
133 
134 
135 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
static const char nl
Definition: Ostream.H:262
Bound the given scalar field if it has gone unbounded.
bool LTS
Definition: createRDeltaT.H:1
messageStream Info
Execute application functionObjects to post-process existing results.
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.