mirror of
https://github.com/OpenFOAM/OpenFOAM-6.git
synced 2025-12-08 06:57:46 +00:00
The old separate incompressible and compressible libraries have been removed. Most of the commonly used RANS and LES models have been upgraded to the new framework but there are a few missing which will be added over the next few days, in particular the realizable k-epsilon model. Some of the less common incompressible RANS models have been introduced into the new library instantiated for incompressible flow only. If they prove to be generally useful they can be templated for compressible and multiphase application. The Spalart-Allmaras DDES and IDDES models have been thoroughly debugged, removing serious errors concerning the use of S rather than Omega. The compressible instances of the models have been augmented by a simple backward-compatible eddyDiffusivity model for thermal transport based on alphat and alphaEff. This will be replaced with a separate run-time selectable thermal transport model framework in a few weeks. For simplicity and ease of maintenance and further development the turbulent transport and wall modeling is based on nut/nuEff rather than mut/muEff for compressible models so that all forms of turbulence models can use the same wall-functions and other BCs. All turbulence model selection made in the constant/turbulenceProperties dictionary with RAS and LES as sub-dictionaries rather than in separate files which added huge complexity for multiphase. All tutorials have been updated so study the changes and update your own cases by comparison with similar cases provided. Sorry for the inconvenience in the break in backward-compatibility but this update to the turbulence modeling is an essential step in the future of OpenFOAM to allow more models to be added and maintained for a wider range of cases and physics. Over the next weeks and months more turbulence models will be added of single and multiphase flow, more additional sub-models and further development and testing of existing models. I hope this brings benefits to all OpenFOAM users. Henry G. Weller
211 lines
6.2 KiB
C
211 lines
6.2 KiB
C
/*---------------------------------------------------------------------------*\
|
|
========= |
|
|
\\ / 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
|
|
PDRFoam
|
|
|
|
Description
|
|
Solver for compressible premixed/partially-premixed combustion with
|
|
turbulence modelling.
|
|
|
|
Combusting RANS code using the b-Xi two-equation model.
|
|
Xi may be obtained by either the solution of the Xi transport
|
|
equation or from an algebraic exression. Both approaches are
|
|
based on Gulder's flame speed correlation which has been shown
|
|
to be appropriate by comparison with the results from the
|
|
spectral model.
|
|
|
|
Strain effects are incorporated directly into the Xi equation
|
|
but not in the algebraic approximation. Further work need to be
|
|
done on this issue, particularly regarding the enhanced removal rate
|
|
caused by flame compression. Analysis using results of the spectral
|
|
model will be required.
|
|
|
|
For cases involving very lean Propane flames or other flames which are
|
|
very strain-sensitive, a transport equation for the laminar flame
|
|
speed is present. This equation is derived using heuristic arguments
|
|
involving the strain time scale and the strain-rate at extinction.
|
|
the transport velocity is the same as that for the Xi equation.
|
|
|
|
For large flames e.g. explosions additional modelling for the flame
|
|
wrinkling due to surface instabilities may be applied.
|
|
|
|
PDR (porosity/distributed resistance) modelling is included to handle
|
|
regions containing blockages which cannot be resolved by the mesh.
|
|
|
|
\*---------------------------------------------------------------------------*/
|
|
|
|
#include "fvCFD.H"
|
|
#include "dynamicFvMesh.H"
|
|
#include "psiuReactionThermo.H"
|
|
#include "turbulentFluidThermoModel.H"
|
|
#include "laminarFlameSpeed.H"
|
|
#include "XiModel.H"
|
|
#include "PDRDragModel.H"
|
|
#include "ignition.H"
|
|
#include "Switch.H"
|
|
#include "bound.H"
|
|
#include "dynamicRefineFvMesh.H"
|
|
#include "pimpleControl.H"
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
int main(int argc, char *argv[])
|
|
{
|
|
#include "setRootCase.H"
|
|
|
|
#include "createTime.H"
|
|
#include "createDynamicFvMesh.H"
|
|
#include "readCombustionProperties.H"
|
|
#include "readGravitationalAcceleration.H"
|
|
#include "createFields.H"
|
|
#include "initContinuityErrs.H"
|
|
#include "readTimeControls.H"
|
|
#include "compressibleCourantNo.H"
|
|
#include "setInitialDeltaT.H"
|
|
|
|
pimpleControl pimple(mesh);
|
|
|
|
scalar StCoNum = 0.0;
|
|
|
|
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
|
|
|
|
Info<< "\nStarting time loop\n" << endl;
|
|
|
|
bool hasChanged = false;
|
|
|
|
while (runTime.run())
|
|
{
|
|
#include "readTimeControls.H"
|
|
#include "compressibleCourantNo.H"
|
|
#include "setDeltaT.H"
|
|
|
|
// Indicators for refinement. Note: before runTime++
|
|
// only for postprocessing reasons.
|
|
tmp<volScalarField> tmagGradP = mag(fvc::grad(p));
|
|
volScalarField normalisedGradP
|
|
(
|
|
"normalisedGradP",
|
|
tmagGradP()/max(tmagGradP())
|
|
);
|
|
normalisedGradP.writeOpt() = IOobject::AUTO_WRITE;
|
|
tmagGradP.clear();
|
|
|
|
runTime++;
|
|
|
|
Info<< "\n\nTime = " << runTime.timeName() << endl;
|
|
|
|
{
|
|
// Make the fluxes absolute
|
|
fvc::makeAbsolute(phi, rho, U);
|
|
|
|
// Test : disable refinement for some cells
|
|
PackedBoolList& protectedCell =
|
|
refCast<dynamicRefineFvMesh>(mesh).protectedCell();
|
|
|
|
if (protectedCell.empty())
|
|
{
|
|
protectedCell.setSize(mesh.nCells());
|
|
protectedCell = 0;
|
|
}
|
|
|
|
forAll(betav, cellI)
|
|
{
|
|
if (betav[cellI] < 0.99)
|
|
{
|
|
protectedCell[cellI] = 1;
|
|
}
|
|
}
|
|
|
|
// Flux estimate for introduced faces.
|
|
volVectorField rhoU("rhoU", rho*U);
|
|
|
|
// Do any mesh changes
|
|
bool meshChanged = mesh.update();
|
|
|
|
|
|
if (meshChanged)
|
|
{
|
|
hasChanged = true;
|
|
}
|
|
|
|
if (runTime.write() && hasChanged)
|
|
{
|
|
betav.write();
|
|
Lobs.write();
|
|
CT.write();
|
|
drag->writeFields();
|
|
flameWrinkling->writeFields();
|
|
hasChanged = false;
|
|
}
|
|
|
|
// Make the fluxes relative to the mesh motion
|
|
fvc::makeRelative(phi, rho, U);
|
|
}
|
|
|
|
|
|
#include "rhoEqn.H"
|
|
|
|
// --- Pressure-velocity PIMPLE corrector loop
|
|
while (pimple.loop())
|
|
{
|
|
#include "UEqn.H"
|
|
|
|
|
|
// --- Pressure corrector loop
|
|
while (pimple.correct())
|
|
{
|
|
#include "bEqn.H"
|
|
#include "ftEqn.H"
|
|
#include "huEqn.H"
|
|
#include "hEqn.H"
|
|
|
|
if (!ign.ignited())
|
|
{
|
|
hu == h;
|
|
}
|
|
|
|
#include "pEqn.H"
|
|
}
|
|
|
|
if (pimple.turbCorr())
|
|
{
|
|
turbulence->correct();
|
|
}
|
|
}
|
|
|
|
runTime.write();
|
|
|
|
Info<< "\nExecutionTime = "
|
|
<< runTime.elapsedCpuTime()
|
|
<< " s\n" << endl;
|
|
}
|
|
|
|
Info<< "\n end\n";
|
|
|
|
return 0;
|
|
}
|
|
|
|
|
|
// ************************************************************************* //
|