openfoam/applications/solvers/lagrangian/porousExplicitSourceReactingParcelFoam/createFields.H
Henry 63da3e9afc Thermodynamics: Rationalization
At the specie level:
    hs = sensible enthalpy
    ha = absolute (what was total) enthalpy
    es = sensibly internal energy
    ea = absolute (what was total) internal energy

At top-level
    Rename total enthalpy h -> ha
    Rename sensible enthalpy hs -> h

Combined h, hs, e and es thermo packages into a single structure.

Thermo packages now provide "he" function which may return either enthalpy or
internal energy, sensible or absolute according to the run-time selected form

alphaEff now returns the effective diffusivity for the particular energy which
the thermodynamics package is selected to solve for.
2012-05-30 15:19:38 +01:00

107 lines
2.5 KiB
C

Info<< "Creating combustion model\n" << endl;
autoPtr<combustionModels::rhoChemistryCombustionModel> combustion
(
combustionModels::rhoChemistryCombustionModel::New
(
mesh
)
);
rhoChemistryModel& chemistry = combustion->pChemistry();
rhoReactionThermo& thermo = chemistry.thermo();
SLGThermo slgThermo(mesh, thermo);
basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();
const word inertSpecie(thermo.lookup("inertSpecie"));
if (!composition.contains(inertSpecie))
{
FatalErrorIn(args.executable())
<< "Specified inert specie '" << inertSpecie << "' not found in "
<< "species list. Available species:" << composition.species()
<< exit(FatalError);
}
volScalarField& p = thermo.p();
volScalarField& hs = thermo.he();
const volScalarField& T = thermo.T();
const volScalarField& psi = thermo.psi();
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
thermo.rho()
);
Info<< "\nReading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
#include "compressibleCreatePhi.H"
Info<< "Creating turbulence model\n" << endl;
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
phi,
thermo
)
);
// Set the turbulence into the combustion model
combustion->setTurbulence(turbulence());
Info<< "Creating multi-variate interpolation scheme\n" << endl;
multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
forAll(Y, i)
{
fields.add(Y[i]);
}
fields.add(hs);
volScalarField dQ
(
IOobject
(
"dQ",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("dQ", dimEnergy/dimTime, 0.0)
);
Info<< "Creating field dpdt\n" << endl;
volScalarField dpdt("dpdt", fvc::ddt(p));
Info<< "Creating field kinetic energy K\n" << endl;
volScalarField K("K", 0.5*magSqr(U));