solidAbsorption.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) 2015-2018 OpenCFD Ltd.
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 "solidAbsorption.H"
30#include "mappedPatchBase.H"
31#include "radiationModel.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38 namespace radiation
39 {
41
43 (
47 );
48 }
49}
50
51
52// * * * * * * * * * * * * * * * * Private members * * * * * * * * * * * * * //
53
54const Foam::fvMesh& Foam::radiation::solidAbsorption::nbrRegion() const
55{
56 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
57 return (refCast<const fvMesh>(mpp.sampleMesh()));
58}
59
60
61Foam::label Foam::radiation::solidAbsorption::nbrPatchIndex() const
62{
63 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
64 return (mpp.samplePolyPatch().index());
65}
66
67
68// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69
71(
72 const dictionary& dict,
73 const polyPatch& pp
74)
75:
77{
78 if (!isA<mappedPatchBase>(pp))
79 {
81 << "\n patch type '" << pp.type()
82 << "' not type '" << mappedPatchBase::typeName << "'"
83 << "\n for patch " << pp.name()
85 }
86}
87
88
89// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
90
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
98(
99 const label bandI,
100 vectorField* incomingDirection,
102) const
103{
104 // Since we're inside initEvaluate/evaluate there might be processor
105 // comms underway. Change the tag we use.
106 int oldTag = UPstream::msgType();
107 UPstream::msgType() = oldTag+1;
108
109 const fvMesh& nbrMesh = nbrRegion();
110
113 (
114 "radiationProperties"
115 );
116
117 scalarField absorptivity
118 (
119 radiation.absorptionEmission().a
120 (
121 bandI
122 )().boundaryField()
123 [
124 nbrPatchIndex()
125 ]
126 );
127
128 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
129
130 mpp.distribute(absorptivity);
131
132 // Restore tag
133 UPstream::msgType() = oldTag;
134
135 return tmp<scalarField>::New(std::move(absorptivity));
136}
137
138
140(
141 const label faceI,
142 const label bandI,
143 const vector dir,
144 const scalar T
145) const
146{
147 return a(bandI, nullptr, nullptr)()[faceI];
148}
149
151(
152 const label bandI,
153 vectorField* incomingDirection,
155) const
156{
157
158 // Since we're inside initEvaluate/evaluate there might be processor
159 // comms underway. Change the tag we use.
160 int oldTag = UPstream::msgType();
161 UPstream::msgType() = oldTag+1;
162
163 const fvMesh& nbrMesh = nbrRegion();
164
167 (
168 "radiationProperties"
169 );
170
171 scalarField emissivity
172 (
173 radiation.absorptionEmission().e
174 (
175 bandI
176 )().boundaryField()
177 [
178 nbrPatchIndex()
179 ]
180 );
181
182 const mappedPatchBase& mpp = refCast<const mappedPatchBase>(pp_);
183
184 mpp.distribute(emissivity);
185
186 // Restore tag
187 UPstream::msgType() = oldTag;
188
189 return tmp<scalarField>::New(std::move(emissivity));
190}
191
192
194(
195 const label faceI,
196 const label bandI,
197 const vector dir,
198 const scalar T
199) const
200{
201 return e(bandI, nullptr, nullptr)()[faceI];
202}
203
204
206{
207 const fvMesh& nbrMesh = nbrRegion();
208
211 (
212 "radiationProperties"
213 );
214
215 return (radiation.absorptionEmission().nBands());
216}
217// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const STLpoint & a() const
Definition: STLtriangleI.H:68
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const Type & lookupObject(const word &name, const bool recursive=false) const
const word & name() const noexcept
The patch name.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
Top level model for radiation modelling.
Radiation absorptivity-emissivity model to be used on walls on inter-region patches when the solid op...
virtual ~solidAbsorption()
Destructor.
label nBands() const
Number of bands.
Based class for wall absorption emission models.
const polyPatch & pp_
Reference to the polyPatch.
A class for managing temporary objects.
Definition: tmp.H:65
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
const dimensionedScalar e
Elementary charge.
Namespace for OpenFOAM.
errorManip< error > abort(error &err)
Definition: errorManip.H:144
IOerror FatalIOError
dictionary dict
volScalarField & e
Definition: createFields.H:11
static const char *const typeName
The type name used in ensight case files.