This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/combustion/dieselEngineFoam/createFields.H

112 lines
2.1 KiB
C++
Raw Permalink Normal View History

Info<< nl << "Reading thermophysicalProperties" << endl;
2010-08-26 14:22:03 +00:00
autoPtr<psiChemistryModel> pChemistry
(
2010-08-26 14:22:03 +00:00
psiChemistryModel::New(mesh)
);
2010-08-26 14:22:03 +00:00
psiChemistryModel& chemistry = pChemistry();
hsCombustionThermo& thermo = chemistry.thermo();
2010-08-26 14:22:03 +00:00
basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();
2010-08-26 14:22:03 +00:00
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 rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh
),
2010-08-26 14:22:03 +00:00
thermo.rho()
);
Info<< "Reading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
2010-08-26 14:22:03 +00:00
volScalarField& p = thermo.p();
const volScalarField& psi = thermo.psi();
const volScalarField& T = thermo.T();
volScalarField& hs = thermo.hs();
#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;
2010-08-26 14:22:03 +00:00
autoPtr<compressible::turbulenceModel> turbulence
(
2010-08-26 14:22:03 +00:00
compressible::turbulenceModel::New
(
rho,
U,
phi,
2010-08-26 14:22:03 +00:00
thermo
)
);
Info<< "Creating field DpDt\n" << endl;
2010-08-26 14:22:03 +00:00
volScalarField DpDt =
fvc::DDt(surfaceScalarField("phiU", phi/fvc::interpolate(rho)), p);
multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
2010-08-26 14:22:03 +00:00
forAll(Y, i)
{
fields.add(Y[i]);
}
2010-08-26 14:22:03 +00:00
fields.add(hs);
2010-08-26 14:22:03 +00:00
DimensionedField<scalar, volMesh> chemistrySh
(
IOobject
(
2010-08-26 14:22:03 +00:00
"chemistry::Sh",
runTime.timeName(),
mesh,
IOobject::NO_READ,
2010-08-26 14:22:03 +00:00
IOobject::NO_WRITE
),
mesh,
2010-08-26 14:22:03 +00:00
dimensionedScalar("chemistrySh", dimEnergy/dimTime/dimVolume, 0.0)
);
mesh.schemesDict().setFluxRequired(p.name());