primaryRadiation.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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2012-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "primaryRadiation.H"
29#include "volFields.H"
31
32// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33
34namespace Foam
35{
36namespace regionModels
37{
38namespace surfaceFilmModels
39{
40
41// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42
44
46(
50);
51
52// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53
55(
57 const dictionary& dict
58)
59:
60 filmRadiationModel(typeName, film, dict),
61 qinPrimary_
62 (
64 (
65 "qin", // same name as qin on primary region to enable mapping
66 film.time().timeName(),
67 film.regionMesh(),
68 IOobject::NO_READ,
69 IOobject::NO_WRITE
70 ),
71 film.regionMesh(),
73 film.mappedPushedFieldPatchTypes<scalar>()
74 )
75{}
76
77
78// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79
81{}
82
83
84// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
85
87{
88 // Transfer qin from primary region
89 qinPrimary_.correctBoundaryConditions();
90}
91
92
94{
96 (
98 (
100 (
101 typeName + ":Shs",
102 film().time().timeName(),
103 film().regionMesh(),
106 ),
107 film().regionMesh(),
109 )
110 );
111
112 scalarField& Shs = tShs.ref();
113 const scalarField& qinP = qinPrimary_;
115
116 Shs = qinP*alpha;
117
118 return tShs;
119}
120
121
122// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
123
124} // End namespace surfaceFilmModels
125} // End namespace regionModels
126} // End namespace Foam
127
128// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
Radiation model whereby the radiative heat flux is mapped from the primary region.
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
volScalarField & alpha
dictionary dict