standardRadiation.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) 2011-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 "standardRadiation.H"
29#include "volFields.H"
32
33// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace regionModels
38{
39namespace surfaceFilmModels
40{
41
42// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43
45
47(
51);
52
53// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54
56(
58 const dictionary& dict
59)
60:
61 filmRadiationModel(typeName, film, dict),
62 qinPrimary_
63 (
65 (
66 "qin", // same name as qin on primary region to enable mapping
67 film.time().timeName(),
68 film.regionMesh(),
69 IOobject::NO_READ,
70 IOobject::NO_WRITE
71 ),
72 film.regionMesh(),
74 film.mappedPushedFieldPatchTypes<scalar>()
75 ),
76 qrNet_
77 (
79 (
80 "qrNet",
81 film.time().timeName(),
82 film.regionMesh(),
83 IOobject::NO_READ,
84 IOobject::NO_WRITE
85 ),
86 film.regionMesh(),
88 zeroGradientFvPatchScalarField::typeName
89 ),
90 beta_(coeffDict_.get<scalar>("beta")),
91 kappaBar_(coeffDict_.get<scalar>("kappaBar"))
92{}
93
94
95// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
96
98{
99 // Transfer qr from primary region
100 qinPrimary_.correctBoundaryConditions();
101}
102
103
105{
107 (
109 (
111 (
112 typeName + ":Shs",
113 film().time().timeName(),
114 film().regionMesh(),
117 ),
118 film().regionMesh(),
120 )
121 );
122
123 scalarField& Shs = tShs.ref();
124 const scalarField& qinP = qinPrimary_;
127
128 Shs = beta_*qinP*alpha*(1.0 - exp(-kappaBar_*delta));
129
130 // Update net qr on local region
131 qrNet_.primitiveFieldRef() = qinP - Shs;
133
134 return tShs;
135}
136
137
138// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
139
140} // End namespace surfaceFilmModels
141} // End namespace regionModels
142} // End namespace Foam
143
144// ************************************************************************* //
scalar delta
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
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.
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
virtual const volScalarField & delta() const =0
Return the film thickness [m].
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 exp(const dimensionedScalar &ds)
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