coldEngineFoam.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-2021 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  Solver for cold-flow in internal combustion engines.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "fvCFD.H"
33 #include "engineTime.H"
34 #include "engineMesh.H"
35 #include "psiThermo.H"
38 #include "OFstream.H"
39 #include "fvModels.H"
40 #include "fvConstraints.H"
41 #include "pimpleControl.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #define CREATE_TIME createEngineTime.H
48  #define CREATE_MESH createEngineMesh.H
49  #include "postProcess.H"
50 
51  #include "setRootCaseLists.H"
52  #include "createEngineTime.H"
53  #include "createEngineMesh.H"
54  #include "createControl.H"
55  #include "createFields.H"
56  #include "createFieldRefs.H"
57  #include "createRhoUf.H"
58  #include "initContinuityErrs.H"
60  #include "compressibleCourantNo.H"
61  #include "setInitialDeltaT.H"
62  #include "startSummary.H"
63 
64  turbulence->validate();
65 
66  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
67 
68  Info<< "\nStarting time loop\n" << endl;
69 
70  while (pimple.run(runTime))
71  {
72  #include "readEngineTimeControls.H"
73  #include "compressibleCourantNo.H"
74  #include "setDeltaT.H"
75 
76  runTime++;
77 
78  Info<< "Engine time = " << runTime.theta() << runTime.unit()
79  << endl;
80 
82 
83  mesh.move();
84 
85  #include "rhoEqn.H"
86 
87  // --- Pressure-velocity PIMPLE corrector loop
88  while (pimple.loop())
89  {
90  fvModels.correct();
91 
92  #include "UEqn.H"
93 
94  // --- Pressure corrector loop
95  while (pimple.correct())
96  {
97  #include "EEqn.H"
98  #include "pEqn.H"
99  }
100 
101  if (pimple.turbCorr())
102  {
103  turbulence->correct();
104  thermophysicalTransport->correct();
105  }
106  }
107 
108  runTime.write();
109 
110  #include "logSummary.H"
111 
112  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
113  << " ClockTime = " << runTime.elapsedClockTime() << " s"
114  << nl << endl;
115  }
116 
117  Info<< "End\n" << endl;
118 
119  return 0;
120 }
121 
122 
123 // ************************************************************************* //
pimpleNoLoopControl & pimple
virtual void correct()
Correct the fvModels.
Definition: fvModels.C:316
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void preUpdateMesh()
Prepare for mesh update.
Definition: fvModels.C:240
dynamicFvMesh & mesh
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
fluidReactionThermophysicalTransportModel & thermophysicalTransport
Foam::fvModels & fvModels
messageStream Info
Execute application functionObjects to post-process existing results.
Creates and initialises the velocity velocity field rhoUf.