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-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  buoyantPimpleFoam
26 
27 Description
28  Transient solver for buoyant, turbulent flow of compressible fluids for
29  ventilation and heat-transfer.
30 
31  Turbulence is modelled using a run-time selectable compressible RAS or
32  LES model.
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "rhoThermo.H"
39 #include "radiationModel.H"
40 #include "fvOptions.H"
41 #include "pimpleControl.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCase.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createFields.H"
54  #include "createFieldRefs.H"
55  #include "createFvOptions.H"
56  #include "initContinuityErrs.H"
57  #include "createTimeControls.H"
58  #include "compressibleCourantNo.H"
59  #include "setInitialDeltaT.H"
60 
61  turbulence->validate();
62 
63  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
64 
65  Info<< "\nStarting time loop\n" << endl;
66 
67  while (runTime.run())
68  {
69  #include "readTimeControls.H"
70  #include "compressibleCourantNo.H"
71  #include "setDeltaT.H"
72 
73  runTime++;
74 
75  Info<< "Time = " << runTime.timeName() << nl << endl;
76 
77  #include "rhoEqn.H"
78 
79  // --- Pressure-velocity PIMPLE corrector loop
80  while (pimple.loop())
81  {
82  #include "UEqn.H"
83  #include "EEqn.H"
84 
85  // --- Pressure corrector loop
86  while (pimple.correct())
87  {
88  #include "pEqn.H"
89  }
90 
91  if (pimple.turbCorr())
92  {
93  turbulence->correct();
94  }
95  }
96 
97  rho = thermo.rho();
98 
99  runTime.write();
100 
101  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
102  << " ClockTime = " << runTime.elapsedClockTime() << " s"
103  << nl << endl;
104  }
105 
106  Info<< "End\n" << endl;
107 
108  return 0;
109 }
110 
111 
112 // ************************************************************************* //
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
psiReactionThermo & thermo
Definition: createFields.H:31
static const char nl
Definition: Ostream.H:262
Info<< "Reading field U\"<< 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\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.