From 42b3f1c9dc79270362639efac9eb7f09f3dc38d7 Mon Sep 17 00:00:00 2001 From: Henry Weller <http://cfd.direct> Date: Fri, 30 Oct 2015 14:32:26 +0000 Subject: [PATCH] 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 --- .../wall/wallHeatFlux/createFields.H | 16 ++++++ .../wall/wallHeatFlux/wallHeatFlux.C | 52 ++++++++++++++++--- 2 files changed, 60 insertions(+), 8 deletions(-) diff --git a/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H b/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H index da7585b8317..e88cc7b37f8 100644 --- a/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H +++ b/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H @@ -42,6 +42,7 @@ if (isA<fluidThermo>(thermo())) const volVectorField& U = UPtr(); #include "compressibleCreatePhi.H" + // Copy phi to autoPtr. Rename to make sure copy is now registered as 'phi'. phi.rename("phiFluid"); phiPtr.reset(new surfaceScalarField("phi", phi)); @@ -54,3 +55,18 @@ if (isA<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) +); diff --git a/applications/utilities/postProcessing/wall/wallHeatFlux/wallHeatFlux.C b/applications/utilities/postProcessing/wall/wallHeatFlux/wallHeatFlux.C index 8c80f0cf217..76550b96993 100644 --- a/applications/utilities/postProcessing/wall/wallHeatFlux/wallHeatFlux.C +++ b/applications/utilities/postProcessing/wall/wallHeatFlux/wallHeatFlux.C @@ -70,19 +70,24 @@ int main(int argc, char *argv[]) const surfaceScalarField::GeometricBoundaryField& patchHeatFlux = heatFlux.boundaryField(); + const volScalarField::GeometricBoundaryField& patchRadHeatFlux = + Qr.boundaryField(); + + const surfaceScalarField::GeometricBoundaryField& magSf = + mesh.magSf().boundaryField(); + Info<< "\nWall heat fluxes [W]" << endl; forAll(patchHeatFlux, patchi) { if (isA<wallFvPatch>(mesh.boundary()[patchi])) { - Info<< mesh.boundary()[patchi].name() - << " " - << gSum - ( - mesh.magSf().boundaryField()[patchi] - *patchHeatFlux[patchi] - ) - << endl; + scalar convFlux = gSum(magSf[patchi]*patchHeatFlux[patchi]); + scalar radFlux = gSum(magSf[patchi]*patchRadHeatFlux[patchi]); + + Info<< mesh.boundary()[patchi].name() << endl + << " convective: " << convFlux << endl + << " radiative: " << radFlux << endl + << " total: " << convFlux + radFlux << endl; } } Info<< endl; @@ -105,6 +110,36 @@ int main(int argc, char *argv[]) } 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; @@ -112,4 +147,5 @@ int main(int argc, char *argv[]) return 0; } + // ************************************************************************* // -- GitLab