From 602c5875e2c9bff97209ebda5a72ade9f009e248 Mon Sep 17 00:00:00 2001 From: andy <andy> Date: Fri, 21 Oct 2011 11:37:57 +0100 Subject: [PATCH] BUG: chemFoam: Calculate the specific gas constant at each time step - Mantis bug: 315 --- .../combustion/chemFoam/createFields.H | 20 +++++++++++++++++++ .../solvers/combustion/chemFoam/pEqn.H | 12 +++++++++-- .../chemFoam/readInitialConditions.H | 3 ++- 3 files changed, 32 insertions(+), 3 deletions(-) diff --git a/applications/solvers/combustion/chemFoam/createFields.H b/applications/solvers/combustion/chemFoam/createFields.H index 4a971492c9c..b37a8818e31 100644 --- a/applications/solvers/combustion/chemFoam/createFields.H +++ b/applications/solvers/combustion/chemFoam/createFields.H @@ -45,9 +45,29 @@ ), thermo.rho() ); + volScalarField& p = thermo.p(); volScalarField& hs = thermo.hs(); + volScalarField Rspecific + ( + IOobject + ( + "Rspecific", + runTime.timeName(), + runTime, + IOobject::NO_READ, + IOobject::AUTO_WRITE + ), + mesh, + dimensionedScalar + ( + "zero", + dimensionSet(dimEnergy/dimMass/dimTemperature), + 0.0 + ) + ); + volVectorField U ( IOobject diff --git a/applications/solvers/combustion/chemFoam/pEqn.H b/applications/solvers/combustion/chemFoam/pEqn.H index 28d240940b6..13f3d603ae7 100644 --- a/applications/solvers/combustion/chemFoam/pEqn.H +++ b/applications/solvers/combustion/chemFoam/pEqn.H @@ -3,7 +3,15 @@ rho = thermo.rho(); if (constProp == "volume") { - p[0] = rho0*R0*thermo.T()[0]; + scalar invW = 0.0; + forAll(Y, i) + { + invW += Y[i][0]/specieData[i].W(); + } + + Rspecific[0] = 1000.0*constant::physicoChemical::R.value()*invW; + + p[0] = rho0*Rspecific[0]*thermo.T()[0]; rho[0] = rho0; } -} \ No newline at end of file +} diff --git a/applications/solvers/combustion/chemFoam/readInitialConditions.H b/applications/solvers/combustion/chemFoam/readInitialConditions.H index 6b26e5fc37e..817ed264f9a 100644 --- a/applications/solvers/combustion/chemFoam/readInitialConditions.H +++ b/applications/solvers/combustion/chemFoam/readInitialConditions.H @@ -106,6 +106,7 @@ scalar rho0 = rho[0]; scalar u0 = hs0 - p0/rho0; scalar R0 = p0/(rho0*T0); - + Rspecific[0] = R0; + scalar integratedHeat = 0.0; -- GitLab