mirror of
https://develop.openfoam.com/Development/openfoam.git
synced 2025-12-28 03:37:59 +00:00
Combined 'dQ()' and 'Sh()' into 'Qdot()' which returns the heat-release rate in the normal units [kg/m/s3] and used as the heat release rate source term in the energy equations, to set the field 'Qdot' in several combustion solvers and for the evaluation of the local time-step when running LTS.
169 lines
3.1 KiB
C
169 lines
3.1 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)
|
|
);
|
|
|
|
|
|
Info<< "Creating field dpdt\n" << endl;
|
|
volScalarField dpdt
|
|
(
|
|
IOobject
|
|
(
|
|
"dpdt",
|
|
runTime.timeName(),
|
|
mesh
|
|
),
|
|
mesh,
|
|
dimensionedScalar("dpdt", p.dimensions()/dimTime, 0)
|
|
);
|
|
|
|
Info<< "Creating field kinetic energy K\n" << endl;
|
|
volScalarField K("K", 0.5*magSqr(U));
|
|
|
|
#include "createClouds.H"
|
|
#include "createSurfaceFilmModel.H"
|
|
#include "createPyrolysisModel.H"
|
|
#include "createRadiationModel.H"
|