wallHeatFlux.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  wallHeatFlux
26 
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.
31 
32 \*---------------------------------------------------------------------------*/
33 
34 #include "fvCFD.H"
36 #include "solidThermo.H"
37 #include "wallFvPatch.H"
38 
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 
41 int main(int argc, char *argv[])
42 {
43  timeSelector::addOptions();
44  #include "addRegionOption.H"
45  #include "setRootCase.H"
46  #include "createTime.H"
47  instantList timeDirs = timeSelector::select0(runTime, args);
48  #include "createNamedMesh.H"
49 
50  forAll(timeDirs, timeI)
51  {
52  runTime.setTime(timeDirs[timeI], timeI);
53  Info<< "Time = " << runTime.timeName() << endl;
54  mesh.readUpdate();
55 
56  #include "createFields.H"
57 
58  surfaceScalarField heatFlux
59  (
61  (
62  (
63  turbulence.valid()
64  ? turbulence->alphaEff()()
65  : thermo->alpha()
66  )
67  )*fvc::snGrad(h)
68  );
69 
70  const surfaceScalarField::Boundary& patchHeatFlux =
71  heatFlux.boundaryField();
72 
73  const volScalarField::Boundary& patchRadHeatFlux =
74  Qr.boundaryField();
75 
76  const surfaceScalarField::Boundary& magSf =
77  mesh.magSf().boundaryField();
78 
79  Info<< "\nWall heat fluxes [W]" << endl;
80  forAll(patchHeatFlux, patchi)
81  {
82  if (isA<wallFvPatch>(mesh.boundary()[patchi]))
83  {
84  scalar convFlux = gSum(magSf[patchi]*patchHeatFlux[patchi]);
85  scalar radFlux = -gSum(magSf[patchi]*patchRadHeatFlux[patchi]);
86 
87  Info<< mesh.boundary()[patchi].name() << endl
88  << " convective: " << convFlux << endl
89  << " radiative: " << radFlux << endl
90  << " total: " << convFlux + radFlux << endl;
91  }
92  }
93  Info<< endl;
94 
95  volScalarField wallHeatFlux
96  (
97  IOobject
98  (
99  "wallHeatFlux",
100  runTime.timeName(),
101  mesh
102  ),
103  mesh,
104  dimensionedScalar("wallHeatFlux", heatFlux.dimensions(), 0.0)
105  );
106 
107  volScalarField::Boundary& wallHeatFluxBf =
108  wallHeatFlux.boundaryFieldRef();
109 
110  forAll(wallHeatFluxBf, patchi)
111  {
112  wallHeatFluxBf[patchi] = patchHeatFlux[patchi];
113  }
114 
115  wallHeatFlux.write();
116 
117  // Write the total heat-flux including the radiative contribution
118  // if available
119  if (Qr.headerOk())
120  {
121  volScalarField totalWallHeatFlux
122  (
123  IOobject
124  (
125  "totalWallHeatFlux",
126  runTime.timeName(),
127  mesh
128  ),
129  mesh,
131  (
132  "totalWallHeatFlux",
133  heatFlux.dimensions(),
134  0.0
135  )
136  );
137 
138  volScalarField::Boundary& totalWallHeatFluxBf =
139  totalWallHeatFlux.boundaryFieldRef();
140 
141  forAll(totalWallHeatFluxBf, patchi)
142  {
143  totalWallHeatFluxBf[patchi] =
144  patchHeatFlux[patchi] - patchRadHeatFlux[patchi];
145  }
146 
147  totalWallHeatFlux.write();
148  }
149  }
150 
151  Info<< "End" << endl;
152 
153  return 0;
154 }
155 
156 
157 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
List< instant > instantList
List of instants.
Definition: instantList.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< surfaceScalarField > interpolate(const RhoType &rho)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dynamicFvMesh & mesh
Type gSum(const FieldField< Field, Type > &f)
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const dimensionedScalar h
Planck constant.
Definition: createFields.H:6
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
messageStream Info
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
volScalarField Qr(IOobject("Qr", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("Qr", dimMass/pow3(dimTime), 0.0))
Foam::argList args(argc, argv)
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45