diff --git a/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H b/applications/utilities/postProcessing/wall/wallHeatFlux/createFields.H
index da7585b83175febac3f111af927b95a21266040e..e88cc7b37f84979b2f08c0471c71f5ef04a9efd2 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 8c80f0cf2177d95c9605846d9cb043d7f7a271da..76550b96993253bf95ec7ca6a9c8e9fcef7b838e 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;
 }
 
+
 // ************************************************************************* //