Files
openfoam/applications/solvers/combustion/fireFoam/fireFoam.C
Mark Olesen 5d9e278e92 ENH: consolidate handling of mandatory/optional command arguments
- for some special cases we wish to mark command-line arguments as
  being optional, in order to do our own treatment. For example,
  when an arbitrary number of arguments should be allowed.

  Now tag this situation with argList::noMandatoryArgs().
  The argList::argsMandatory() query can then be used in any further
  logic, including the standard default argument checking.

- with the new default check, can consolidate the special-purpose

      "setRootCaseNonMandatoryArgs.H"

  into the regular

      "setRootCase.H"

- revert to a simple "setRootCase.H" and move all the listing related
  bits to a "setRootCaseLists.H" file. This leaves the information
  available for solvers, or whoever else wishes, without being
  introduced everywhere.

- add include guards and scoping to the listing files and rename to
  something less generic.

     listOptions.H -> setRootCaseListOptions.H
     listOutput.H  -> setRootCaseListOutput.H
2018-12-13 01:45:09 +01:00

138 lines
3.7 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2017 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Application
fireFoam
Group
grpCombustionSolvers
Description
Transient solver for fires and turbulent diffusion flames with reacting
particle clouds, surface film and pyrolysis modelling.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "turbulentFluidThermoModel.H"
#include "basicReactingCloud.H"
#include "surfaceFilmModel.H"
#include "pyrolysisModelCollection.H"
#include "radiationModel.H"
#include "SLGThermo.H"
#include "solidChemistryModel.H"
#include "psiReactionThermo.H"
#include "CombustionModel.H"
#include "pimpleControl.H"
#include "fvOptions.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
argList::addNote
(
"Transient solver for fires and turbulent diffusion flames"
" with reacting particle clouds, surface film and pyrolysis modelling."
);
#include "postProcess.H"
#include "addCheckCaseOptions.H"
#include "setRootCaseLists.H"
#include "createTime.H"
#include "createMesh.H"
#include "createControl.H"
#include "createFields.H"
#include "createFieldRefs.H"
#include "initContinuityErrs.H"
#include "createTimeControls.H"
#include "compressibleCourantNo.H"
#include "setInitialDeltaT.H"
#include "createRegionControls.H"
turbulence->validate();
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info<< "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readTimeControls.H"
#include "compressibleCourantNo.H"
#include "solidRegionDiffusionNo.H"
#include "setMultiRegionDeltaT.H"
#include "setDeltaT.H"
++runTime;
Info<< "Time = " << runTime.timeName() << nl << endl;
parcels.evolve();
surfaceFilm.evolve();
if(solvePyrolysisRegion)
{
pyrolysis.evolve();
}
if (solvePrimaryRegion)
{
#include "rhoEqn.H"
// --- PIMPLE loop
while (pimple.loop())
{
#include "UEqn.H"
#include "YEEqn.H"
// --- Pressure corrector loop
while (pimple.correct())
{
#include "pEqn.H"
}
if (pimple.turbCorr())
{
turbulence->correct();
}
}
rho = thermo.rho();
}
runTime.write();
runTime.printExecutionTime(Info);
}
Info<< "End" << endl;
return 0;
}
// ************************************************************************* //