Loading...
Searching...
No Matches
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 Copyright (C) 2023 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "standardRadiation.H"
30#include "volFields.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace regionModels
39{
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
46
48(
49 filmRadiationModel,
50 standardRadiation,
52);
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
56standardRadiation::standardRadiation
57(
59 const dictionary& dict
60)
61:
62 filmRadiationModel(typeName, film, dict),
63 qinPrimary_
64 (
66 (
67 "qin", // same name as qin on primary region to enable mapping
68 film.time().timeName(),
69 film.regionMesh().thisDb(),
70 IOobject::NO_READ,
71 IOobject::NO_WRITE,
72 IOobject::REGISTER
73 ),
74 film.regionMesh(),
76 film.mappedPushedFieldPatchTypes<scalar>()
77 ),
78 qrNet_
79 (
81 (
82 "qrNet",
83 film.time().timeName(),
84 film.regionMesh().thisDb(),
85 IOobject::NO_READ,
86 IOobject::NO_WRITE,
87 IOobject::REGISTER
88 ),
89 film.regionMesh(),
91 fvPatchFieldBase::zeroGradientType()
92 ),
93 beta_(coeffDict_.get<scalar>("beta")),
94 kappaBar_(coeffDict_.get<scalar>("kappaBar"))
95{}
96
97
98// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
99
101{
102 // Transfer qr from primary region
103 qinPrimary_.correctBoundaryConditions();
104}
105
106
108{
109 auto tShs = volScalarField::New
110 (
113 film().regionMesh(),
115 );
116 scalarField& Shs = tShs.ref();
117
118 const scalarField& qinP = qinPrimary_;
121
122 Shs = beta_*qinP*alpha*(1.0 - exp(-kappaBar_*delta));
123
124 // Update net qr on local region
125 qrNet_.primitiveFieldRef() = qinP - Shs;
127
128 return tShs;
129}
130
131
132// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133
134} // End namespace surfaceFilmModels
135} // End namespace regionModels
136} // End namespace Foam
137
138// ************************************************************************* //
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
void correctBoundaryConditions()
Correct boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Template invariant parts for fvPatchField.
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].
const dictionary coeffDict_
Coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
A class for managing temporary objects.
Definition tmp.H:75
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
word timeName
Definition getTimeIndex.H:3
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
volScalarField & alpha