buoyantPimpleFoam.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  buoyantPimpleFoam
26 
27 Group
28  grpHeatTransferSolvers
29 
30 Description
31  Transient solver for buoyant, turbulent flow of compressible fluids for
32  ventilation and heat-transfer.
33 
34  Turbulence is modelled using a run-time selectable compressible RAS or
35  LES model.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "fvCFD.H"
40 #include "rhoThermo.H"
42 #include "radiationModel.H"
43 #include "fvOptions.H"
44 #include "pimpleControl.H"
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
51  #include "setRootCase.H"
52  #include "createTime.H"
53  #include "createMesh.H"
54 
55  pimpleControl pimple(mesh);
56 
57  #include "createFields.H"
58  #include "createMRF.H"
59  #include "createFvOptions.H"
60  #include "createRadiationModel.H"
61  #include "initContinuityErrs.H"
62  #include "createTimeControls.H"
63  #include "compressibleCourantNo.H"
64  #include "setInitialDeltaT.H"
65 
66  turbulence->validate();
67 
68  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70  Info<< "\nStarting time loop\n" << endl;
71 
72  while (runTime.run())
73  {
74  #include "createTimeControls.H"
75  #include "compressibleCourantNo.H"
76  #include "setDeltaT.H"
77 
78  runTime++;
79 
80  Info<< "Time = " << runTime.timeName() << nl << endl;
81 
82  #include "rhoEqn.H"
83 
84  // --- Pressure-velocity PIMPLE corrector loop
85  while (pimple.loop())
86  {
87  #include "UEqn.H"
88  #include "EEqn.H"
89 
90  // --- Pressure corrector loop
91  while (pimple.correct())
92  {
93  #include "pEqn.H"
94  }
95 
96  if (pimple.turbCorr())
97  {
98  turbulence->correct();
99  }
100  }
101 
102  rho = thermo.rho();
103 
104  runTime.write();
105 
106  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
107  << " ClockTime = " << runTime.elapsedClockTime() << " s"
108  << nl << endl;
109  }
110 
111  Info<< "End\n" << endl;
112 
113  return 0;
114 }
115 
116 
117 // ************************************************************************* //
fvOptions.H
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
createRadiationModel.H
pimpleControl.H
createFvOptions.H
rhoThermo.H
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
rho
rho
Definition: pEqn.H:3
setRootCase.H
createTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
createTime.H
fvCFD.H
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
fixedFluxPressureFvPatchScalarField.H
radiationModel.H
turbulentFluidThermoModel.H