initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / wall / wallHeatFlux / wallHeatFlux.C
blobc83fcb983ac18b8280536b1f4f6e136a98a66666
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software; you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     wallHeatFlux
28 Description
29     Calculates and writes the heat flux for all patches as the boundary field
30     of a volScalarField and also prints the integrated flux for all wall
31     patches.
33 \*---------------------------------------------------------------------------*/
35 #include "fvCFD.H"
36 #include "hCombustionThermo.H"
37 #include "RASModel.H"
38 #include "wallFvPatch.H"
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 int main(int argc, char *argv[])
44     timeSelector::addOptions();
45     #include "setRootCase.H"
46     #include "createTime.H"
47     instantList timeDirs = timeSelector::select0(runTime, args);
48     #include "createMesh.H"
50     forAll(timeDirs, timeI)
51     {
52         runTime.setTime(timeDirs[timeI], timeI);
53         Info<< "Time = " << runTime.timeName() << endl;
54         mesh.readUpdate();
56         #include "createFields.H"
58         surfaceScalarField heatFlux =
59             fvc::interpolate(RASModel->alphaEff())*fvc::snGrad(h);
61         const surfaceScalarField::GeometricBoundaryField& patchHeatFlux =
62             heatFlux.boundaryField();
64         Info<< "\nWall heat fluxes [W]" << endl;
65         forAll(patchHeatFlux, patchi)
66         {
67             if (typeid(mesh.boundary()[patchi]) == typeid(wallFvPatch))
68             {
69                 Info<< mesh.boundary()[patchi].name()
70                     << " "
71                     << sum
72                        (
73                            mesh.magSf().boundaryField()[patchi]
74                           *patchHeatFlux[patchi]
75                        )
76                     << endl;
77             }
78         }
79         Info<< endl;
81         volScalarField wallHeatFlux
82         (
83             IOobject
84             (
85                 "wallHeatFlux",
86                 runTime.timeName(),
87                 mesh
88             ),
89             mesh,
90             dimensionedScalar("wallHeatFlux", heatFlux.dimensions(), 0.0)
91         );
93         forAll(wallHeatFlux.boundaryField(), patchi)
94         {
95             wallHeatFlux.boundaryField()[patchi] = patchHeatFlux[patchi];
96         }
98         wallHeatFlux.write();
99     }
101     Info<< "End" << endl;
103     return 0;
106 // ************************************************************************* //