standardRadiation.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 \*---------------------------------------------------------------------------*/
25 
26 #include "standardRadiation.H"
27 #include "volFields.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace regionModels
36 {
37 namespace surfaceFilmModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(standardRadiation, 0);
43 
45 (
46  filmRadiationModel,
47  standardRadiation,
48  dictionary
49 );
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
56  const dictionary& dict
57 )
58 :
59  filmRadiationModel(typeName, film, dict),
60  qinPrimary_
61  (
62  IOobject
63  (
64  "qin", // same name as qin on primary region to enable mapping
65  film.time().timeName(),
66  film.regionMesh(),
69  ),
70  film.regionMesh(),
72  film.mappedPushedFieldPatchTypes<scalar>()
73  ),
74  qrNet_
75  (
76  IOobject
77  (
78  "qrNet",
79  film.time().timeName(),
80  film.regionMesh(),
83  ),
84  film.regionMesh(),
86  zeroGradientFvPatchScalarField::typeName
87  ),
88  beta_(readScalar(coeffDict_.lookup("beta"))),
89  kappaBar_(readScalar(coeffDict_.lookup("kappaBar")))
90 {}
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
100 
102 {
103  // Transfer qr from primary region
104  qinPrimary_.correctBoundaryConditions();
105 }
106 
107 
109 {
111  (
113  (
114  typeName + ":Shs",
115  film().regionMesh(),
117  )
118  );
119 
120  scalarField& Shs = tShs.ref();
121  const scalarField& qinP = qinPrimary_;
122  const scalarField& delta = filmModel_.delta();
123  const scalarField& alpha = filmModel_.alpha();
124 
125  Shs = beta_*qinP*alpha*(1.0 - exp(-kappaBar_*delta));
126 
127  // Update net qr on local region
128  qrNet_.primitiveFieldRef() = qinP - Shs;
129  qrNet_.correctBoundaryConditions();
130 
131  return tShs;
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
136 
137 } // End namespace surfaceFilmModels
138 } // End namespace regionModels
139 } // End namespace Foam
140 
141 // ************************************************************************* //
scalar delta
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
const Time & time() const
Return the reference to the time database.
Definition: regionModelI.H:37
addToRunTimeSelectionTable(surfaceFilmRegionModel, kinematicSingleLayer, mesh)
Macros for easy insertion into run-time selection tables.
const surfaceFilmRegionModel & film() const
Return const access to the film surface film model.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
dimensionedScalar exp(const dimensionedScalar &ds)
virtual const volScalarField & delta() const =0
Return the film thickness [m].
const fvMesh & regionMesh() const
Return the region mesh database.
Definition: regionModelI.H:61
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if successful.
Definition: doubleScalar.H:68
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
dimensionedScalar pow3(const dimensionedScalar &ds)
wordList mappedPushedFieldPatchTypes() const
Return boundary types for pushed mapped field patches.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
void correctBoundaryConditions()
Correct boundary field.
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Namespace for OpenFOAM.
standardRadiation(surfaceFilmRegionModel &film, const dictionary &dict)
Construct from surface film model and dictionary.