Files
openfoam/applications/solvers/combustion/fireFoam/createFields.H
2017-06-21 12:20:58 +01:00

156 lines
2.8 KiB
C

Info<< "Creating combustion model\n" << endl;
autoPtr<combustionModels::psiCombustionModel> combustion
(
combustionModels::psiCombustionModel::New
(
mesh
)
);
Info<< "Reading thermophysical properties\n" << endl;
psiReactionThermo& thermo = combustion->thermo();
thermo.validate(args.executable(), "h", "e");
SLGThermo slgThermo(mesh, thermo);
basicMultiComponentMixture& composition = thermo.composition();
PtrList<volScalarField>& Y = composition.Y();
const word inertSpecie(thermo.lookup("inertSpecie"));
if (!composition.species().found(inertSpecie))
{
FatalIOErrorIn(args.executable().c_str(), thermo)
<< "Inert specie " << inertSpecie << " not found in available species "
<< composition.species()
<< exit(FatalIOError);
}
Info<< "Creating field rho\n" << endl;
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
thermo.rho()
);
volScalarField& p = thermo.p();
Info<< "\nReading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
#include "compressibleCreatePhi.H"
#include "createMRF.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());
#include "readGravitationalAcceleration.H"
#include "readhRef.H"
#include "gh.H"
#include "readpRef.H"
volScalarField p_rgh
(
IOobject
(
"p_rgh",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
mesh.setFluxRequired(p_rgh.name());
#include "phrghEqn.H"
multivariateSurfaceInterpolationScheme<scalar>::fieldTable fields;
forAll(Y, i)
{
fields.add(Y[i]);
}
fields.add(thermo.he());
IOdictionary additionalControlsDict
(
IOobject
(
"additionalControls",
runTime.constant(),
mesh,
IOobject::MUST_READ_IF_MODIFIED,
IOobject::NO_WRITE
)
);
Switch solvePrimaryRegion
(
additionalControlsDict.lookup("solvePrimaryRegion")
);
Switch solvePyrolysisRegion
(
additionalControlsDict.lookupOrDefault<bool>("solvePyrolysisRegion", true)
);
volScalarField Qdot
(
IOobject
(
"Qdot",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
mesh,
dimensionedScalar("Qdot", dimEnergy/dimVolume/dimTime, 0.0)
);
#include "createDpdt.H"
#include "createK.H"
#include "createClouds.H"
#include "createSurfaceFilmModel.H"
#include "createPyrolysisModel.H"
#include "createRadiationModel.H"