wallHeatFlux: Add support for radiative and total heat-fluxes

Patch provided by Daniel Jasinski
Resolves feature request http://www.openfoam.org/mantisbt/view.php?id=1856
This commit is contained in:
Henry Weller
2015-10-30 14:32:26 +00:00
parent 941117a1c4
commit 5deddf7e30
2 changed files with 60 additions and 8 deletions

View File

@ -42,6 +42,7 @@ if (isA<fluidThermo>(thermo()))
const volVectorField& U = UPtr(); const volVectorField& U = UPtr();
#include "compressibleCreatePhi.H" #include "compressibleCreatePhi.H"
// Copy phi to autoPtr. Rename to make sure copy is now registered as 'phi'. // Copy phi to autoPtr. Rename to make sure copy is now registered as 'phi'.
phi.rename("phiFluid"); phi.rename("phiFluid");
phiPtr.reset(new surfaceScalarField("phi", phi)); phiPtr.reset(new surfaceScalarField("phi", phi));
@ -54,3 +55,18 @@ if (isA<fluidThermo>(thermo()))
refCast<const fluidThermo>(thermo()) refCast<const fluidThermo>(thermo())
); );
} }
// Read radiative heat-flux if available
volScalarField Qr
(
IOobject
(
"Qr",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::NO_WRITE
),
mesh,
dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0)
);

View File

@ -70,19 +70,24 @@ int main(int argc, char *argv[])
const surfaceScalarField::GeometricBoundaryField& patchHeatFlux = const surfaceScalarField::GeometricBoundaryField& patchHeatFlux =
heatFlux.boundaryField(); heatFlux.boundaryField();
const volScalarField::GeometricBoundaryField& patchRadHeatFlux =
Qr.boundaryField();
const surfaceScalarField::GeometricBoundaryField& magSf =
mesh.magSf().boundaryField();
Info<< "\nWall heat fluxes [W]" << endl; Info<< "\nWall heat fluxes [W]" << endl;
forAll(patchHeatFlux, patchi) forAll(patchHeatFlux, patchi)
{ {
if (isA<wallFvPatch>(mesh.boundary()[patchi])) if (isA<wallFvPatch>(mesh.boundary()[patchi]))
{ {
Info<< mesh.boundary()[patchi].name() scalar convFlux = gSum(magSf[patchi]*patchHeatFlux[patchi]);
<< " " scalar radFlux = gSum(magSf[patchi]*patchRadHeatFlux[patchi]);
<< gSum
( Info<< mesh.boundary()[patchi].name() << endl
mesh.magSf().boundaryField()[patchi] << " convective: " << convFlux << endl
*patchHeatFlux[patchi] << " radiative: " << radFlux << endl
) << " total: " << convFlux + radFlux << endl;
<< endl;
} }
} }
Info<< endl; Info<< endl;
@ -105,6 +110,36 @@ int main(int argc, char *argv[])
} }
wallHeatFlux.write(); wallHeatFlux.write();
// Write the total heat-flux including the radiative contribution
// if available
if (Qr.headerOk())
{
volScalarField totalWallHeatFlux
(
IOobject
(
"totalWallHeatFlux",
runTime.timeName(),
mesh
),
mesh,
dimensionedScalar
(
"totalWallHeatFlux",
heatFlux.dimensions(),
0.0
)
);
forAll(totalWallHeatFlux.boundaryField(), patchi)
{
totalWallHeatFlux.boundaryField()[patchi] =
patchHeatFlux[patchi] + patchRadHeatFlux[patchi];
}
totalWallHeatFlux.write();
}
} }
Info<< "End" << endl; Info<< "End" << endl;
@ -112,4 +147,5 @@ int main(int argc, char *argv[])
return 0; return 0;
} }
// ************************************************************************* // // ************************************************************************* //