BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / wall / wallHeatFlux / wallHeatFlux.C
blob070869ce28e93fa995d493afdc475613af847e75
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     wallHeatFlux
27 Description
28     Calculates and writes the heat flux for all patches as the boundary field
29     of a volScalarField and also prints the integrated flux for all wall
30     patches.
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "hCombustionThermo.H"
36 #include "RASModel.H"
37 #include "wallFvPatch.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43     timeSelector::addOptions();
44     #include "setRootCase.H"
45     #include "createTime.H"
46     instantList timeDirs = timeSelector::select0(runTime, args);
47     #include "createMesh.H"
49     forAll(timeDirs, timeI)
50     {
51         runTime.setTime(timeDirs[timeI], timeI);
52         Info<< "Time = " << runTime.timeName() << endl;
53         mesh.readUpdate();
55         #include "createFields.H"
57         surfaceScalarField heatFlux
58         (
59             fvc::interpolate(RASModel->alphaEff())*fvc::snGrad(h)
60         );
62         const surfaceScalarField::GeometricBoundaryField& patchHeatFlux =
63             heatFlux.boundaryField();
65         Info<< "\nWall heat fluxes [W]" << endl;
66         forAll(patchHeatFlux, patchi)
67         {
68             if (isA<wallFvPatch>(mesh.boundary()[patchi]))
69             {
70                 Info<< mesh.boundary()[patchi].name()
71                     << " "
72                     << sum
73                        (
74                            mesh.magSf().boundaryField()[patchi]
75                           *patchHeatFlux[patchi]
76                        )
77                     << endl;
78             }
79         }
80         Info<< endl;
82         volScalarField wallHeatFlux
83         (
84             IOobject
85             (
86                 "wallHeatFlux",
87                 runTime.timeName(),
88                 mesh
89             ),
90             mesh,
91             dimensionedScalar("wallHeatFlux", heatFlux.dimensions(), 0.0)
92         );
94         forAll(wallHeatFlux.boundaryField(), patchi)
95         {
96             wallHeatFlux.boundaryField()[patchi] = patchHeatFlux[patchi];
97         }
99         wallHeatFlux.write();
100     }
102     Info<< "End" << endl;
104     return 0;
107 // ************************************************************************* //