All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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-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  Transient solver for compressible, turbulent engine flow with a spray
29  particle cloud.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "engineMesh.H"
37 #include "combustionModel.H"
38 #include "fvModels.H"
39 #include "fvConstraints.H"
40 #include "pimpleControl.H"
41 #include "pressureReference.H"
42 #include "CorrectPhi.H"
43 #include "localEulerDdtScheme.H"
44 #include "fvcSmooth.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 "createEngineControls.H"
58  #include "initContinuityErrs.H"
59  #include "createFields.H"
60  #include "createFieldRefs.H"
61  #include "createRhoUfIfPresent.H"
62  #include "startSummary.H"
63 
64  turbulence->validate();
65 
66  if (!LTS)
67  {
68  #include "compressibleCourantNo.H"
69  #include "setInitialDeltaT.H"
70  }
71 
72  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74  Info<< "\nStarting time loop\n" << endl;
75 
76  while (pimple.run(runTime))
77  {
78  #include "readEngineControls.H"
79 
80  // Store divrhoU from the previous mesh so that it can be mapped
81  // and used in correctPhi to ensure the corrected phi has the
82  // same divergence
83  autoPtr<volScalarField> divrhoU;
84  if (correctPhi)
85  {
86  divrhoU = new volScalarField
87  (
88  "divrhoU",
90  );
91  }
92 
93  if (LTS)
94  {
95  #include "setRDeltaT.H"
96  }
97  else
98  {
99  #include "compressibleCourantNo.H"
100  #include "setDeltaT.H"
101  }
102 
103  runTime++;
104 
105  Info<< "Engine time = " << runTime.theta() << runTime.unit() << endl;
106 
107  // Store momentum to set rhoUf for introduced faces.
108  autoPtr<volVectorField> rhoU;
109  if (rhoUf.valid())
110  {
111  rhoU = new volVectorField("rhoU", rho*U);
112  }
113 
115 
116  // Do any mesh changes
117  mesh.move();
118 
119  if (mesh.changing())
120  {
121  MRF.update();
122 
123  if (correctPhi)
124  {
125  #include "../../../compressible/rhoPimpleFoam/correctPhi.H"
126  }
127 
128  if (checkMeshCourantNo)
129  {
130  #include "meshCourantNo.H"
131  }
132  }
133 
134  if (!pimple.simpleRho())
135  {
136  #include "rhoEqn.H"
137  }
138 
139  // --- PIMPLE loop
140  while (pimple.loop())
141  {
142  fvModels.correct();
143 
144  #include "UEqn.H"
145  #include "YEqn.H"
146  #include "EEqn.H"
147 
148  // --- Pressure corrector loop
149  while (pimple.correct())
150  {
151  #include "../../../compressible/rhoPimpleFoam/pEqn.H"
152  }
153 
154  if (pimple.turbCorr())
155  {
156  turbulence->correct();
157  thermophysicalTransport->correct();
158  }
159  }
160 
161  #include "logSummary.H"
162 
163  rho = thermo.rho();
164 
165  runTime.write();
166 
167  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
168  << " ClockTime = " << runTime.elapsedClockTime() << " s"
169  << nl << endl;
170  }
171 
172  Info<< "End" << endl;
173 
174  return 0;
175 }
176 
177 
178 // ************************************************************************* //
fluidReactionThermo & thermo
Definition: createFields.H:28
pimpleNoLoopControl & pimple
IOMRFZoneList & MRF
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
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
correctPhi
checkMeshCourantNo
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
Creates and initialises the velocity field rhoUf if required.
dynamicFvMesh & mesh
phi
Definition: correctPhi.H:3
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
bool LTS
Definition: createRDeltaT.H:1
fluidReactionThermophysicalTransportModel & thermophysicalTransport
Foam::fvModels & fvModels
Calculates and outputs the mean and maximum Courant Numbers.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
messageStream Info
Execute application functionObjects to post-process existing results.
autoPtr< surfaceVectorField > rhoUf
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...