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