sprayFoam.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-2018 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  sprayFoam
26 
27 Description
28  Transient solver for compressible, turbulent flow with a spray particle
29  cloud.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
35 #include "basicSprayCloud.H"
36 #include "psiReactionThermo.H"
37 #include "CombustionModel.H"
38 #include "radiationModel.H"
39 #include "SLGThermo.H"
40 #include "pimpleControl.H"
41 #include "fvOptions.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 int main(int argc, char *argv[])
46 {
47  #include "postProcess.H"
48 
49  #include "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createControl.H"
53  #include "createTimeControls.H"
54  #include "createFields.H"
55  #include "createFieldRefs.H"
56  #include "compressibleCourantNo.H"
57  #include "setInitialDeltaT.H"
58  #include "initContinuityErrs.H"
59 
60  turbulence->validate();
61 
62  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
63 
64  Info<< "\nStarting time loop\n" << endl;
65 
66  while (runTime.run())
67  {
68  #include "readTimeControls.H"
69  #include "compressibleCourantNo.H"
70  #include "setDeltaT.H"
71 
72  runTime++;
73 
74  Info<< "Time = " << runTime.timeName() << nl << endl;
75 
76  parcels.evolve();
77 
78  if (!pimple.frozenFlow())
79  {
80  #include "rhoEqn.H"
81 
82  // --- Pressure-velocity PIMPLE corrector loop
83  while (pimple.loop())
84  {
85  #include "UEqn.H"
86  #include "YEqn.H"
87  #include "EEqn.H"
88 
89  // --- Pressure corrector loop
90  while (pimple.correct())
91  {
92  #include "pEqn.H"
93  }
94 
95  if (pimple.turbCorr())
96  {
97  turbulence->correct();
98  }
99  }
100 
101  rho = thermo.rho();
102 
103  if (runTime.write())
104  {
105  combustion->Qdot()().write();
106  }
107  }
108  else
109  {
110  if (runTime.writeTime())
111  {
112  parcels.write();
113  }
114  }
115 
116  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
117  << " ClockTime = " << runTime.elapsedClockTime() << " s"
118  << nl << endl;
119  }
120 
121  Info<< "End\n" << endl;
122 
123  return 0;
124 }
125 
126 
127 // ************************************************************************* //
pimpleNoLoopControl & pimple
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Read the control parameters used by setDeltaT.
rhoReactionThermo & thermo
Definition: createFields.H:28
void write(Ostream &, const label, const dictionary &)
Write with dictionary lookup.
static const char nl
Definition: Ostream.H:265
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
Info<< "Creating combustion model\"<< endl;autoPtr< CombustionModel< psiReactionThermo > > combustion(CombustionModel< psiReactionThermo >::New(thermo, turbulence()))
messageStream Info
Execute application functionObjects to post-process existing results.
Read the control parameters used by setDeltaT.