This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/utilities/postProcessing/wall/wallHeatFlux/wallHeatFlux.C
2016-06-21 15:04:12 +02:00

104 lines
3.2 KiB
C

/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | foam-extend: Open Source CFD
\\ / O peration | Version: 4.0
\\ / A nd | Web: http://www.foam-extend.org
\\/ M anipulation | For copyright notice see file Copyright
-------------------------------------------------------------------------------
License
This file is part of foam-extend.
foam-extend is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by the
Free Software Foundation, either version 3 of the License, or (at your
option) any later version.
foam-extend is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
General Public License for more details.
You should have received a copy of the GNU General Public License
along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
Application
wallHeatFlux
Description
Calculates and writes the heat flux for all patches as the boundary field
of a volScalarField and also prints the integrated flux for all wall
patches.
\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "hCombustionThermo.H"
#include "RASModel.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
{
timeSelector::addOptions();
#include "setRootCase.H"
#include "createTime.H"
instantList timeDirs = timeSelector::select0(runTime, args);
#include "createMesh.H"
forAll(timeDirs, timeI)
{
runTime.setTime(timeDirs[timeI], timeI);
Info<< "Time = " << runTime.timeName() << endl;
mesh.readUpdate();
#include "createFields.H"
surfaceScalarField heatFlux =
fvc::interpolate(RASModel->alphaEff())*fvc::snGrad(h);
const surfaceScalarField::GeometricBoundaryField& patchHeatFlux =
heatFlux.boundaryField();
Info<< "\nWall heat fluxes [W]" << endl;
forAll(patchHeatFlux, patchi)
{
if (mesh.boundary()[patchi].isWall())
{
Info<< mesh.boundary()[patchi].name()
<< " "
<< sum
(
mesh.magSf().boundaryField()[patchi]
*patchHeatFlux[patchi]
)
<< endl;
}
}
Info<< endl;
volScalarField wallHeatFlux
(
IOobject
(
"wallHeatFlux",
runTime.timeName(),
mesh
),
mesh,
dimensionedScalar("wallHeatFlux", heatFlux.dimensions(), 0.0)
);
forAll(wallHeatFlux.boundaryField(), patchi)
{
wallHeatFlux.boundaryField()[patchi] = patchHeatFlux[patchi];
}
wallHeatFlux.write();
}
Info<< "End" << endl;
return 0;
}
// ************************************************************************* //