From 03adb37ab4ab988cf04ae72ab5fcdb508905eef0 Mon Sep 17 00:00:00 2001
From: Henry Weller <http://cfd.direct>
Date: Mon, 15 Jun 2015 22:36:27 +0100
Subject: [PATCH] LTSCoalChemistryFoam: LTS version of coalChemistryFoam

---
 .../LTSCoalChemistryFoam.C                    | 125 ++++++++++++++++++
 .../LTSCoalChemistryFoam/Make/files           |   3 +
 .../LTSCoalChemistryFoam/Make/options         |  58 ++++++++
 .../LTSCoalChemistryFoam/timeScales.H         | 114 ++++++++++++++++
 .../lagrangian/coalChemistryFoam/Make/options |   2 -
 5 files changed, 300 insertions(+), 2 deletions(-)
 create mode 100644 applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/LTSCoalChemistryFoam.C
 create mode 100755 applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/files
 create mode 100755 applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/options
 create mode 100644 applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/timeScales.H

diff --git a/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/LTSCoalChemistryFoam.C b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/LTSCoalChemistryFoam.C
new file mode 100644
index 00000000000..5ee21e9c49a
--- /dev/null
+++ b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/LTSCoalChemistryFoam.C
@@ -0,0 +1,125 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
+   \\    /   O peration     |
+    \\  /    A nd           | Copyright (C) 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
+    LTSCoalChemistryFoam
+
+Description
+    Local time stepping (LTS) solver for steady simulation of:
+    - 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 "fvcSmooth.H"
+
+// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+int main(int argc, char *argv[])
+{
+    #include "setRootCase.H"
+
+    #include "createTime.H"
+    #include "createMesh.H"
+
+    pimpleControl pimple(mesh);
+
+    #include "readGravitationalAcceleration.H"
+    #include "createFields.H"
+    #include "readTimeControls.H"
+    #include "createRDeltaT.H"
+    #include "createMRF.H"
+    #include "createFvOptions.H"
+    #include "createClouds.H"
+    #include "createRadiationModel.H"
+    #include "initContinuityErrs.H"
+
+    // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
+
+    Info<< "\nStarting time loop\n" << endl;
+
+    while (runTime.run())
+    {
+        #include "readTimeControls.H"
+
+        runTime++;
+
+        Info<< "Time = " << runTime.timeName() << nl << endl;
+
+        rhoEffLagrangian = coalParcels.rhoEff() + limestoneParcels.rhoEff();
+        pDyn = 0.5*rho*magSqr(U);
+
+        coalParcels.evolve();
+
+        limestoneParcels.evolve();
+
+        #include "timeScales.H"
+
+        #include "rhoEqn.H"
+
+        // --- Pressure-velocity PIMPLE corrector loop
+        while (pimple.loop())
+        {
+            if (pimple.turbCorr())
+            {
+                turbulence->correct();
+            }
+
+            #include "UEqn.H"
+            #include "YEqn.H"
+            #include "EEqn.H"
+
+            // --- Pressure corrector loop
+            while (pimple.correct())
+            {
+                #include "pEqn.H"
+            }
+        }
+
+        runTime.write();
+
+        Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
+            << "  ClockTime = " << runTime.elapsedClockTime() << " s"
+            << nl << endl;
+    }
+
+    Info<< "End\n" << endl;
+
+    return(0);
+}
+
+
+// ************************************************************************* //
diff --git a/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/files b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/files
new file mode 100755
index 00000000000..b0d27fc8717
--- /dev/null
+++ b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/files
@@ -0,0 +1,3 @@
+LTSCoalChemistryFoam.C
+
+EXE = $(FOAM_APPBIN)/LTSCoalChemistryFoam
diff --git a/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/options b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/options
new file mode 100755
index 00000000000..d777874d6c1
--- /dev/null
+++ b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/Make/options
@@ -0,0 +1,58 @@
+EXE_INC = \
+    -I.. \
+    -I../../reactingParcelFoam/LTSReactingParcelFoam \
+    -I$(LIB_SRC)/finiteVolume/lnInclude \
+    -I${LIB_SRC}/meshTools/lnInclude \
+    -I$(LIB_SRC)/TurbulenceModels/turbulenceModels/lnInclude \
+    -I$(LIB_SRC)/TurbulenceModels/compressible/lnInclude \
+    -I$(LIB_SRC)/lagrangian/basic/lnInclude \
+    -I$(LIB_SRC)/lagrangian/intermediate/lnInclude \
+    -I$(LIB_SRC)/lagrangian/coalCombustion/lnInclude \
+    -I$(LIB_SRC)/lagrangian/distributionModels/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/specie/lnInclude \
+    -I$(LIB_SRC)/transportModels/compressible/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/properties/liquidProperties/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/properties/liquidMixtureProperties/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/properties/solidProperties/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/properties/solidMixtureProperties/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/thermophysicalFunctions/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/reactionThermo/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/SLGThermo/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/chemistryModel/lnInclude \
+    -I$(LIB_SRC)/thermophysicalModels/radiation/lnInclude \
+    -I$(LIB_SRC)/regionModels/regionModel/lnInclude \
+    -I$(LIB_SRC)/regionModels/surfaceFilmModels/lnInclude \
+    -I$(LIB_SRC)/ODE/lnInclude \
+    -I$(LIB_SRC)/combustionModels/lnInclude \
+    -I$(FOAM_SOLVERS)/combustion/reactingFoam \
+    -I$(LIB_SRC)/fvOptions/lnInclude \
+    -I$(LIB_SRC)/sampling/lnInclude
+
+EXE_LIBS = \
+    -lfiniteVolume \
+    -lmeshTools \
+    -lturbulenceModels \
+    -lcompressibleTurbulenceModels \
+    -llagrangian \
+    -llagrangianIntermediate \
+    -llagrangianTurbulence \
+    -lcoalCombustion\
+    -lspecie \
+    -lcompressibleTransportModels \
+    -lfluidThermophysicalModels \
+    -lliquidProperties \
+    -lliquidMixtureProperties \
+    -lsolidProperties \
+    -lsolidMixtureProperties \
+    -lthermophysicalFunctions \
+    -lreactionThermophysicalModels \
+    -lSLGThermo \
+    -lchemistryModel \
+    -lradiationModels \
+    -lregionModels \
+    -lsurfaceFilmModels \
+    -lODE \
+    -lcombustionModels \
+    -lfvOptions \
+    -lsampling
diff --git a/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/timeScales.H b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/timeScales.H
new file mode 100644
index 00000000000..a32048fc04f
--- /dev/null
+++ b/applications/solvers/lagrangian/coalChemistryFoam/LTSCoalChemistryFoam/timeScales.H
@@ -0,0 +1,114 @@
+/*---------------------------------------------------------------------------*\
+  =========                 |
+  \\      /  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/>.
+
+\*---------------------------------------------------------------------------*/
+
+Info<< "Time scales min/max:" << endl;
+
+{
+    // Cache old time scale field
+    tmp<volScalarField> trDeltaT
+    (
+        new volScalarField
+        (
+            IOobject
+            (
+                "rDeltaT0",
+                runTime.timeName(),
+                mesh,
+                IOobject::NO_READ,
+                IOobject::NO_WRITE,
+                false
+            ),
+            rDeltaT
+        )
+    );
+    const volScalarField& rDeltaT0 = trDeltaT();
+
+
+    // Flow time scale
+    // ~~~~~~~~~~~~~~~
+    {
+        rDeltaT =
+            fvc::surfaceSum
+            (
+                mag(phi)*mesh.deltaCoeffs()/(maxCo*mesh.magSf())
+            )
+           /rho;
+
+        rDeltaT.max(1.0/maxDeltaT);
+
+        Info<< "    Flow        = "
+            << gMin(1/rDeltaT.internalField()) << ", "
+            << gMax(1/rDeltaT.internalField()) << endl;
+    }
+
+
+    // Temperature source time scale
+    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+
+    {
+        scalarField tau
+        (
+            runTime.deltaTValue()
+           *mag
+            (
+                (coalParcels.hsTrans() + limestoneParcels.hsTrans())
+                    /(mesh.V()*runTime.deltaT())
+              + combustion->Sh()()
+              - fvc::div(phi, 0.5*magSqr(U), "div(phi,K)")()
+            )
+           /rho
+        );
+
+        tau = alphaTemp*thermo.Cp()*T/(tau + ROOTVSMALL);
+
+        Info<< "    Temperature = " << min(maxDeltaT, gMin(tau)) << ", "
+            << min(maxDeltaT, gMax(tau)) << endl;
+
+        rDeltaT.internalField() = max(rDeltaT.internalField(), 1/tau);
+    }
+
+
+    // Limit rate of change of time scale
+    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+    // - reduce as much as required for flow, but limit source contributions
+    const dimensionedScalar deltaTRamp("deltaTRamp", dimless, 1/(1 + 0.2));
+    rDeltaT = max(rDeltaT, rDeltaT0*deltaTRamp);
+
+
+    // Limit the largest time scale
+    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+    rDeltaT.max(1/maxDeltaT);
+
+
+    // Spatially smooth the time scale field
+    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
+    fvc::smooth(rDeltaT, rDeltaTSmoothingCoeff);
+
+    Info<< "    Overall     = " << min(1/rDeltaT).value()
+        << ", " << max(1/rDeltaT).value() << nl << endl;
+}
+
+
+// ************************************************************************* //
diff --git a/applications/solvers/lagrangian/coalChemistryFoam/Make/options b/applications/solvers/lagrangian/coalChemistryFoam/Make/options
index f678b0ce3f6..35e1cb7fb61 100644
--- a/applications/solvers/lagrangian/coalChemistryFoam/Make/options
+++ b/applications/solvers/lagrangian/coalChemistryFoam/Make/options
@@ -27,8 +27,6 @@ EXE_INC = \
     -I$(LIB_SRC)/fvOptions/lnInclude \
     -I$(LIB_SRC)/sampling/lnInclude
 
-
-
 EXE_LIBS = \
     -lfiniteVolume \
     -lmeshTools \
-- 
GitLab