constantRadiation.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  Copyright (C) 2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "constantRadiation.H"
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace regionModels
38 {
39 namespace surfaceFilmModels
40 {
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 defineTypeNameAndDebug(constantRadiation, 0);
45 
47 (
48  filmRadiationModel,
51 );
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
55 constantRadiation::constantRadiation
56 (
58  const dictionary& dict
59 )
60 :
61  filmRadiationModel(typeName, film, dict),
62  qrConst_
63  (
64  IOobject
65  (
66  IOobject::scopedName(typeName, "qrConst"),
67  film.time().timeName(),
68  film.regionMesh().thisDb(),
69  IOobject::MUST_READ,
70  IOobject::AUTO_WRITE
71  ),
72  film.regionMesh()
73  ),
74  mask_
75  (
76  IOobject
77  (
78  IOobject::scopedName(typeName, "mask"),
79  film.time().timeName(),
80  film.regionMesh().thisDb(),
81  IOobject::READ_IF_PRESENT,
82  IOobject::AUTO_WRITE
83  ),
84  film.regionMesh(),
85  dimensionedScalar("one", dimless, 1.0)
86  ),
87  absorptivity_(coeffDict_.get<scalar>("absorptivity")),
88  timeStart_(coeffDict_.get<scalar>("timeStart")),
89  duration_(coeffDict_.get<scalar>("duration"))
90 {
91  mask_ = pos0(mask_);
92 }
93 
94 
95 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
96 
98 {}
99 
101 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
102 
104 {}
105 
106 
108 {
109  auto tShs = volScalarField::New
110  (
111  IOobject::scopedName(typeName, "Shs"),
113  film().regionMesh(),
115  );
116  auto& Shs = tShs.ref().primitiveFieldRef();
117 
118  const scalar time = film().time().value();
119 
120  if ((time >= timeStart_) && (time <= timeStart_ + duration_))
121  {
122  const scalarField& qr = qrConst_;
123  const scalarField& alpha = filmModel_.alpha();
124 
125  Shs = mask_*qr*alpha*absorptivity_;
126  }
127 
128  return tShs;
129 }
130 
131 
132 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 
134 } // End namespace surfaceFilmModels
135 } // End namespace regionModels
136 } // End namespace Foam
137 
138 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const Type & value() const noexcept
Return const reference to value.
dictionary dict
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
virtual const volScalarField & alpha() const =0
Return the film coverage, 1 = covered, 0 = uncovered / [].
Film constant radiation model. The constant radiative flux is specified by the user, and operated over a time interval defined by a start time and duration. In addition, a mask can be applied to shield the film from the radiation.
const dimensionSet dimless
Dimensionless.
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.
word timeName
Definition: getTimeIndex.H:3
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
dimensionedScalar pos0(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > Shs()
Return the radiation sensible enthalpy source.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
surfaceFilmRegionModel & filmModel_
Reference to the film surface film model.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
defineTypeNameAndDebug(kinematicSingleLayer, 0)
Do not request registration (bool: false)
const Time & time() const noexcept
Return the reference to the time database.
Definition: regionModel.H:244
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127