Loading...
Searching...
No Matches
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"
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
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{
57 return (refCast<const fvMesh>(mpp.sampleMesh()));
58}
59
60
61Foam::label Foam::radiation::solidAbsorption::nbrPatchIndex() const
62{
64 return (mpp.samplePolyPatch().index());
65}
66
67
68// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
69
71(
72 const dictionary& dict,
73 const polyPatch& pp
74)
75:
77{
79 {
81 << "\n patch type '" << pp.type()
82 << "' not type '" << mappedPatchBase::typeName << "'"
83 << "\n for patch " << pp.name()
86}
87
88
89// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92{}
93
94
95// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
96
98(
99 const label bandI,
100 const vectorField* incomingDirection,
101 const scalarField* T
102) const
103{
104 // Since we're inside initEvaluate/evaluate there might be processor
105 // comms underway. Change the tag we use.
106 const int oldTag = UPstream::incrMsgType();
107
108 const fvMesh& nbrMesh = nbrRegion();
109
112 (
113 "radiationProperties"
114 );
115
116 scalarField absorptivity
117 (
118 radiation.absorptionEmission().a
119 (
120 bandI
121 )().boundaryField()
122 [
123 nbrPatchIndex()
124 ]
125 );
126
128
129 mpp.distribute(absorptivity);
131 UPstream::msgType(oldTag); // Restore tag
132
133 return tmp<scalarField>::New(std::move(absorptivity));
134}
135
136
138(
139 const label faceI,
140 const label bandI,
141 const vector dir,
142 const scalar T
143) const
144{
145 return a(bandI, nullptr, nullptr)()[faceI];
146}
147
149(
150 const label bandI,
151 const vectorField* incomingDirection,
152 const scalarField* T
153) const
154{
155
156 // Since we're inside initEvaluate/evaluate there might be processor
157 // comms underway. Change the tag we use.
158 const int oldTag = UPstream::incrMsgType();
159
160 const fvMesh& nbrMesh = nbrRegion();
161
164 (
165 "radiationProperties"
166 );
167
168 scalarField emissivity
169 (
170 radiation.absorptionEmission().e
171 (
172 bandI
173 )().boundaryField()
174 [
175 nbrPatchIndex()
176 ]
177 );
178
180
181 mpp.distribute(emissivity);
183 UPstream::msgType(oldTag); // Restore tag
184
185 return tmp<scalarField>::New(std::move(emissivity));
186}
187
188
190(
191 const label faceI,
192 const label bandI,
193 const vector dir,
194 const scalar T
195) const
196{
197 return e(bandI, nullptr, nullptr)()[faceI];
198}
199
200
202{
203 const fvMesh& nbrMesh = nbrRegion();
204
207 (
208 "radiationProperties"
209 );
210
211 return (radiation.absorptionEmission().nBands());
212}
213// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Determines a mapping between patch face centres and mesh cell or face centres and processors they're ...
const polyMesh & sampleMesh() const
Get the region mesh.
const polyPatch & samplePolyPatch() const
Get the patch on the region.
void distribute(List< Type > &lst) const
Wrapper around map/interpolate data distribution.
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
label index() const noexcept
The index of this patch in the boundaryMesh.
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
Top level model for radiation modelling.
Radiation absorptivity-emissivity model to be used on walls on inter-region patches when the solid op...
tmp< scalarField > a(const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
absorptivity coefficient
virtual ~solidAbsorption()
Destructor.
label nBands() const
Number of bands.
tmp< scalarField > e(const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Return emission coefficient.
solidAbsorption(const dictionary &dict, const polyPatch &pp)
Construct from components.
Based class for wall absorption emission models.
wallAbsorptionEmissionModel(const dictionary &dict, const polyPatch &pp)
Construct from components.
const polyPatch & pp_
Reference to the polyPatch.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for radiation modelling.
Namespace for OpenFOAM.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
Vector< scalar > vector
Definition vector.H:57
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
volScalarField & e