FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

wallHeatFlux.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     wallHeatFlux
00026 
00027 Description
00028     Calculates and writes the heat flux for all patches as the boundary field
00029     of a volScalarField and also prints the integrated flux for all wall
00030     patches.
00031 
00032 Usage
00033 
00034     - wallHeatFlux [OPTIONS]
00035 
00036     @param -noZero \n
00037     Ignore timestep 0.
00038 
00039     @param -constant \n
00040     Include the constant directory.
00041 
00042     @param -time <time>\n
00043     Apply only to specific time.
00044 
00045     @param -latestTime \n
00046     Only apply to latest time step.
00047 
00048     @param -case <dir>\n
00049     Case directory.
00050 
00051     @param -parallel \n
00052     Run in parallel.
00053 
00054     @param -help \n
00055     Display help message.
00056 
00057     @param -doc \n
00058     Display Doxygen API documentation page for this application.
00059 
00060     @param -srcDoc \n
00061     Display Doxygen source documentation page for this application.
00062 
00063 \*---------------------------------------------------------------------------*/
00064 
00065 #include <finiteVolume/fvCFD.H>
00066 #include <reactionThermophysicalModels/hCombustionThermo.H>
00067 #include <compressibleRASModels/RASModel.H>
00068 #include <finiteVolume/wallFvPatch.H>
00069 
00070 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00071 
00072 int main(int argc, char *argv[])
00073 {
00074     timeSelector::addOptions();
00075     #include <OpenFOAM/setRootCase.H>
00076     #include <OpenFOAM/createTime.H>
00077     instantList timeDirs = timeSelector::select0(runTime, args);
00078     #include <OpenFOAM/createMesh.H>
00079 
00080     forAll(timeDirs, timeI)
00081     {
00082         runTime.setTime(timeDirs[timeI], timeI);
00083         Info<< "Time = " << runTime.timeName() << endl;
00084         mesh.readUpdate();
00085 
00086         #include "createFields.H"
00087 
00088         surfaceScalarField heatFlux =
00089             fvc::interpolate(RASModel->alphaEff())*fvc::snGrad(h);
00090 
00091         const surfaceScalarField::GeometricBoundaryField& patchHeatFlux =
00092             heatFlux.boundaryField();
00093 
00094         Info<< "\nWall heat fluxes [W]" << endl;
00095         forAll(patchHeatFlux, patchi)
00096         {
00097             if (isA<wallFvPatch>(mesh.boundary()[patchi]))
00098             {
00099                 Info<< mesh.boundary()[patchi].name()
00100                     << " "
00101                     << sum
00102                        (
00103                            mesh.magSf().boundaryField()[patchi]
00104                           *patchHeatFlux[patchi]
00105                        )
00106                     << endl;
00107             }
00108         }
00109         Info<< endl;
00110 
00111         volScalarField wallHeatFlux
00112         (
00113             IOobject
00114             (
00115                 "wallHeatFlux",
00116                 runTime.timeName(),
00117                 mesh
00118             ),
00119             mesh,
00120             dimensionedScalar("wallHeatFlux", heatFlux.dimensions(), 0.0)
00121         );
00122 
00123         forAll(wallHeatFlux.boundaryField(), patchi)
00124         {
00125             wallHeatFlux.boundaryField()[patchi] = patchHeatFlux[patchi];
00126         }
00127 
00128         wallHeatFlux.write();
00129     }
00130 
00131     Info<< "End" << endl;
00132 
00133     return 0;
00134 }
00135 
00136 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines