engineFoam.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-2020 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  engineFoam
26 
27 Description
28  Transient solver for compressible, turbulent engine flow with a spray
29  particle cloud.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "engineTime.H"
35 #include "engineMesh.H"
38 #include "basicSprayCloud.H"
39 #include "psiReactionThermo.H"
40 #include "CombustionModel.H"
41 #include "radiationModel.H"
42 #include "SLGThermo.H"
43 #include "pimpleControl.H"
44 #include "fvOptions.H"
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 int main(int argc, char *argv[])
49 {
50  #define CREATE_TIME createEngineTime.H
51  #define CREATE_MESH createEngineMesh.H
52  #include "postProcess.H"
53 
54  #include "setRootCaseLists.H"
55  #include "createEngineTime.H"
56  #include "createEngineMesh.H"
57  #include "createControl.H"
58  #include "readEngineTimeControls.H"
59  #include "createFields.H"
60  #include "createFieldRefs.H"
61  #include "compressibleCourantNo.H"
62  #include "setInitialDeltaT.H"
63  #include "initContinuityErrs.H"
64  #include "createRhoUfIfPresent.H"
65  #include "startSummary.H"
66 
67  turbulence->validate();
68 
69  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
70 
71  Info<< "\nStarting time loop\n" << endl;
72 
73  while (pimple.run(runTime))
74  {
75  #include "readEngineTimeControls.H"
76  #include "compressibleCourantNo.H"
77  #include "setDeltaT.H"
78 
79  runTime++;
80 
81  Info<< "Engine time = " << runTime.theta() << runTime.unit() << endl;
82 
83  mesh.move();
84 
85  parcels.evolve();
86 
87  #include "rhoEqn.H"
88 
89  // --- Pressure-velocity PIMPLE corrector loop
90  while (pimple.loop())
91  {
92  #include "UEqn.H"
93  #include "YEqn.H"
94  #include "EEqn.H"
95 
96  // --- Pressure corrector loop
97  while (pimple.correct())
98  {
99  #include "pEqn.H"
100  }
101 
102  if (pimple.turbCorr())
103  {
104  turbulence->correct();
105  thermophysicalTransport->correct();
106  }
107  }
108 
109  #include "logSummary.H"
110 
111  rho = thermo.rho();
112 
113  runTime.write();
114 
115  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
116  << " ClockTime = " << runTime.elapsedClockTime() << " s"
117  << nl << endl;
118  }
119 
120  Info<< "End\n" << endl;
121 
122  return 0;
123 }
124 
125 
126 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
rhoReactionThermo & thermo
Definition: createFields.H:28
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
rhoReactionThermophysicalTransportModel & thermophysicalTransport
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
messageStream Info
Execute application functionObjects to post-process existing results.