adding rhoReactingFoam solver

This commit is contained in:
andy 2009-07-06 15:33:55 +01:00
parent cd0a0a696c
commit 516c20c8e7
8 changed files with 392 additions and 0 deletions

View File

@ -0,0 +1,3 @@
rhoReactingFoam.C
EXE = $(FOAM_APPBIN)/rhoReactingFoam

View File

@ -0,0 +1,19 @@
EXE_INC = \
-I../XiFoam \
-I$(LIB_SRC)/turbulenceModels/compressible/turbulenceModel \
-I$(LIB_SRC)/thermophysicalModels/specie/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/reactionThermo/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/basic/lnInclude \
-I$(LIB_SRC)/thermophysicalModels/chemistryModel/lnInclude \
-I$(LIB_SRC)/ODE/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude
EXE_LIBS = \
-lcompressibleRASModels \
-lcompressibleLESModels \
-lreactionThermophysicalModels \
-lspecie \
-lbasicThermophysicalModels \
-lchemistryModel \
-lODE \
-lfiniteVolume

View File

@ -0,0 +1,43 @@
tmp<fv::convectionScheme<scalar> > mvConvection
(
fv::convectionScheme<scalar>::New
(
mesh,
fields,
phi,
mesh.divScheme("div(phi,Yi_h)")
)
);
{
label inertIndex = -1;
volScalarField Yt = 0.0*Y[0];
for (label i=0; i<Y.size(); i++)
{
if (Y[i].name() != inertSpecie)
{
volScalarField& Yi = Y[i];
solve
(
fvm::ddt(rho, Yi)
+ mvConvection->fvmDiv(phi, Yi)
- fvm::laplacian(turbulence->muEff(), Yi)
==
kappa*chemistry.RR(i),
mesh.solver("Yi")
);
Yi.max(0.0);
Yt += Yi;
}
else
{
inertIndex = i;
}
}
Y[inertIndex] = scalar(1) - Yt;
Y[inertIndex].max(0.0);
}

View File

@ -0,0 +1,24 @@
{
Info << "Solving chemistry" << endl;
chemistry.solve
(
runTime.value() - runTime.deltaT().value(),
runTime.deltaT().value()
);
// turbulent time scale
if (turbulentReaction)
{
volScalarField tk =
Cmix*sqrt(turbulence->muEff()/rho/turbulence->epsilon());
volScalarField tc = chemistry.tc();
// Chalmers PaSR model
kappa = (runTime.deltaT() + tc)/(runTime.deltaT() + tc + tk);
}
else
{
kappa = 1.0;
}
}

View File

@ -0,0 +1,85 @@
Info<< nl << "Reading thermophysicalProperties" << endl;
autoPtr<rhoChemistryModel> pChemistry
(
rhoChemistryModel::New(mesh)
);
rhoChemistryModel& chemistry = pChemistry();
hReactionThermo& thermo = chemistry.thermo();
basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();
word inertSpecie(thermo.lookup("inertSpecie"));
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh
),
thermo.rho()
);
Info<< "Reading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
volScalarField& p = thermo.p();
const volScalarField& psi = thermo.psi();
volScalarField& h = thermo.h();
#include "compressibleCreatePhi.H"
volScalarField kappa
(
IOobject
(
"kappa",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("zero", dimless, 0.0)
);
Info << "Creating turbulence model.\n" << nl;
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
phi,
thermo
)
);
Info<< "Creating field DpDt\n" << endl;
volScalarField DpDt =
fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
forAll (Y, i)
{
fields.add(Y[i]);
}
fields.add(h);

View File

@ -0,0 +1,93 @@
{
rho = thermo.rho();
// Thermodynamic density needs to be updated by psi*d(p) after the
// pressure solution - done in 2 parts. Part 1:
thermo.rho() -= psi*p;
volScalarField rUA = 1.0/UEqn.A();
U = rUA*UEqn.H();
if (transonic)
{
surfaceScalarField phiv =
(fvc::interpolate(U) & mesh.Sf())
+ fvc::ddtPhiCorr(rUA, rho, U, phi);
phi = fvc::interpolate(rho)*phiv;
surfaceScalarField phid
(
"phid",
fvc::interpolate(thermo.psi())*phiv
);
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix pEqn
(
fvc::ddt(rho) + fvc::div(phi)
+ correction(fvm::ddt(psi, p) + fvm::div(phid, p))
- fvm::laplacian(rho*rUA, p)
);
if (ocorr == nOuterCorr && corr == nCorr && nonOrth == nNonOrthCorr)
{
pEqn.solve(mesh.solver(p.name() + "Final"));
}
else
{
pEqn.solve();
}
if (nonOrth == nNonOrthCorr)
{
phi += pEqn.flux();
}
}
}
else
{
phi =
fvc::interpolate(rho)
*(
(fvc::interpolate(U) & mesh.Sf())
+ fvc::ddtPhiCorr(rUA, rho, U, phi)
);
for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
{
fvScalarMatrix pEqn
(
fvc::ddt(rho) + psi*correction(fvm::ddt(p))
+ fvc::div(phi)
- fvm::laplacian(rho*rUA, p)
);
if (ocorr == nOuterCorr && corr == nCorr && nonOrth == nNonOrthCorr)
{
pEqn.solve(mesh.solver(p.name() + "Final"));
}
else
{
pEqn.solve();
}
if (nonOrth == nNonOrthCorr)
{
phi += pEqn.flux();
}
}
}
// Second part of thermodynamic density update
thermo.rho() += psi*p;
#include "rhoEqn.H"
#include "compressibleContinuityErrs.H"
U -= rUA*fvc::grad(p);
U.correctBoundaryConditions();
DpDt = fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
}

View File

@ -0,0 +1,23 @@
Info<< "Reading chemistry properties\n" << endl;
IOdictionary chemistryProperties
(
IOobject
(
"chemistryProperties",
runTime.constant(),
mesh,
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
)
);
Switch turbulentReaction(chemistryProperties.lookup("turbulentReaction"));
dimensionedScalar Cmix("Cmix", dimless, 1.0);
if (turbulentReaction)
{
chemistryProperties.lookup("Cmix") >> Cmix;
}

View File

@ -0,0 +1,102 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2009-2009 OpenCFD Ltd.
\\/ 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 2 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, write to the Free Software Foundation,
Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Application
rhoReactingFoam
Description
Chemical reaction code using density based thermodynamics package.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "hReactionThermo.H"
#include "turbulenceModel.H"
#include "rhoChemistryModel.H"
#include "chemistrySolver.H"
#include "multivariateScheme.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
#include "setRootCase.H"
#include "createTime.H"
#include "createMesh.H"
#include "readChemistryProperties.H"
#include "readEnvironmentalProperties.H"
#include "createFields.H"
#include "initContinuityErrs.H"
#include "readTimeControls.H"
#include "compressibleCourantNo.H"
#include "setInitialDeltaT.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
Info << "\nStarting time loop\n" << endl;
while (runTime.run())
{
#include "readTimeControls.H"
#include "readPISOControls.H"
#include "compressibleCourantNo.H"
#include "setDeltaT.H"
runTime++;
Info<< "Time = " << runTime.timeName() << nl << endl;
#include "chemistry.H"
#include "rhoEqn.H"
#include "UEqn.H"
for (label ocorr=1; ocorr <= nOuterCorr; ocorr++)
{
#include "YEqn.H"
#include "hEqn.H"
// --- PISO loop
for (int corr=1; corr<=nCorr; corr++)
{
#include "pEqn.H"
}
}
turbulence->correct();
rho = thermo.rho();
runTime.write();
Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
<< " ClockTime = " << runTime.elapsedClockTime() << " s"
<< nl << endl;
}
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //