/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2011-2015 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 coalChemistryFoam Description Transient solver for: - compressible, - turbulent flow, with - coal and limestone parcel injections, - energy source, and - combustion. \*---------------------------------------------------------------------------*/ #include "fvCFD.H" #include "turbulentFluidThermoModel.H" #include "basicThermoCloud.H" #include "coalCloud.H" #include "psiCombustionModel.H" #include "fvIOoptionList.H" #include "radiationModel.H" #include "SLGThermo.H" #include "pimpleControl.H" #include "localEulerDdtScheme.H" #include "fvcSmooth.H" // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // int main(int argc, char *argv[]) { #include "setRootCase.H" #include "createTime.H" #include "createMesh.H" pimpleControl pimple(mesh); #include "createRDeltaT.H" #include "initContinuityErrs.H" #include "readGravitationalAcceleration.H" #include "createFields.H" #include "createMRF.H" #include "createFvOptions.H" #include "createClouds.H" #include "createRadiationModel.H" if (!LTS) { #include "readTimeControls.H" #include "CourantNo.H" #include "setInitialDeltaT.H" } // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // Info<< "\nStarting time loop\n" << endl; while (runTime.run()) { #include "readTimeControls.H" if (LTS) { #include "setRDeltaT.H" } else { #include "compressibleCourantNo.H" #include "setDeltaT.H" } runTime++; Info<< "Time = " << runTime.timeName() << nl << endl; rhoEffLagrangian = coalParcels.rhoEff() + limestoneParcels.rhoEff(); pDyn = 0.5*rho*magSqr(U); coalParcels.evolve(); limestoneParcels.evolve(); #include "rhoEqn.H" // --- Pressure-velocity PIMPLE corrector loop while (pimple.loop()) { #include "UEqn.H" #include "YEqn.H" #include "EEqn.H" // --- Pressure corrector loop while (pimple.correct()) { #include "pEqn.H" } if (pimple.turbCorr()) { turbulence->correct(); } } rho = thermo.rho(); runTime.write(); Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s" << " ClockTime = " << runTime.elapsedClockTime() << " s" << nl << endl; } Info<< "End\n" << endl; return 0; } // ************************************************************************* //