Loading...
Searching...
No Matches
boundaryRadiationProperties.H
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-2022 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
26Class
27 Foam::radiation::boundaryRadiationProperties
28
29Description
30 Boundary radiation properties holder
31
32SourceFiles
33 boundaryRadiationProperties.C
34
35\*---------------------------------------------------------------------------*/
36
37#ifndef Foam_boundaryRadiationProperties_H
38#define Foam_boundaryRadiationProperties_H
39
40#include "MeshObject.H"
42
43// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44
45namespace Foam
46{
47
48// Forward Declarations
49class fvMesh;
50
51namespace radiation
52{
54/*---------------------------------------------------------------------------*\
55 Class boundaryRadiationProperties Declaration
56\*---------------------------------------------------------------------------*/
57
59:
60 public MeshObject
61 <
62 fvMesh,
63 GeometricMeshObject,
64 boundaryRadiationProperties
65 >
66{
67 // Private Typedefs
68
69 typedef MeshObject
70 <
71 fvMesh,
74 > MeshObject_type;
75
76
77 // Private Data
78
79 //- Per patch the boundaryRadiationProperties
81 radBoundaryPropertiesPtrList_;
82
83 //- Per (face)zone the boundaryRadiationProperties
85 radZonePropertiesPtrList_;
86
87
88public:
89
90 // Declare name of the class and its debug switch
91 TypeName("boundaryRadiationProperties");
92
93
94 // Constructors
95
96 //- Construct given fvMesh
98
99
100 // Member Functions
101
102 //- Return identifiers of face zones activated for boundary radiation
103 const labelList faceZoneIds() const
104 {
105 DynamicList<label> ncZones;
107 forAll(radZonePropertiesPtrList_, i)
108 {
109 if (radZonePropertiesPtrList_.test(i))
110 {
111 ncZones.append(i);
112 }
113 }
114
115 return ncZones;
116 }
117
118 //- Access boundary emissivity on patch
120 (
121 const label patchI,
122 const label bandI = 0,
123 const vectorField* incomingDirection = nullptr,
124 const scalarField* T = nullptr
125 ) const;
126
127 //- Access boundary emissivity on face
128 scalar faceEmissivity
129 (
130 const label patchI,
131 const label faceI,
132 const label bandI = 0,
133 vector incomingDirection = Zero,
134 scalar T = 0
135 ) const;
136
137 //- Access boundary absorptivity on patch
139 (
140 const label patchI,
141 const label bandI = 0,
142 const vectorField* incomingDirection = nullptr,
143 const scalarField* T = nullptr
144 ) const;
145
146 //- Access boundary absorptivity on face
147 scalar faceAbsorptivity
148 (
149 const label patchI,
150 const label faceI,
151 const label bandI = 0,
152 vector incomingDirection = Zero,
153 scalar T = 0
154 ) const;
155
156 //- Access boundary transmissivity on patch
158 (
159 const label patchI,
160 const label bandI = 0,
161 const vectorField* incomingDirection = nullptr,
162 const scalarField* T = nullptr
163 ) const;
164
165 //- Access boundary transmissivity on face
166 scalar faceTransmissivity
167 (
168 const label patchI,
169 const label faceI,
170 const label bandI = 0,
171 vector incomingDirection = Zero,
172 scalar T = 0
173 ) const;
174
175 // Specialisations for faceZones
176
177 //- Access transmissivity on set of (internal) faces. Zone name only
178 // used to lookup the properties in boundaryRadiationProperties
180 (
181 const label zoneI,
182 const labelUList& faceIDs, // internal faces
183 const label bandI = 0,
184 vector incomingDirection = Zero,
185 scalar T = 0
186 ) const;
187
188
189 //- Access boundary diffuse reflectivity on patch
191 (
192 const label patchI,
193 const label bandI = 0,
194 const vectorField* incomingDirection = nullptr,
195 const scalarField* T = nullptr
196 ) const;
197
198 //- Access boundary diffuse reflectivity on face
200 (
201 const label patchI,
202 const label faceI,
203 const label bandI = 0,
204 vector incomingDirection = Zero,
205 scalar T = 0
206 ) const;
207
208 //- Access boundary specular reflectivity on patch
210 (
211 const label patchI,
212 const label bandI = 0,
213 const vectorField* incomingDirection = nullptr,
214 const scalarField* T = nullptr
215 ) const;
216
217 //- Access boundary specular reflectivity on face
219 (
220 const label patchI,
221 const label faceI,
222 const label bandI = 0,
223 vector incomingDirection = Zero,
224 scalar T = 0
225 ) const;
226
227
228 //- Destructor
230};
231
232
233// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
234
235} // End namespace radiation
236} // End namespace Foam
237
238// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239
240
241#endif
242
243// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
GeometricMeshObject(const word &objName, const objectRegistry &obr)
Definition MeshObject.H:307
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
scalar faceTransmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary transmissivity on face.
const labelList faceZoneIds() const
Return identifiers of face zones activated for boundary radiation.
scalar faceSpecReflectivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary specular reflectivity on face.
boundaryRadiationProperties(const fvMesh &mesh)
Construct given fvMesh.
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary transmissivity on patch.
TypeName("boundaryRadiationProperties")
scalar faceDiffReflectivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary diffuse reflectivity on face.
scalar faceAbsorptivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary absorptivity on face.
tmp< scalarField > zoneTransmissivity(const label zoneI, const labelUList &faceIDs, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access transmissivity on set of (internal) faces. Zone name only.
tmp< scalarField > absorptivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary absorptivity on patch.
tmp< scalarField > emissivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary emissivity on patch.
tmp< scalarField > specReflectivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary specular reflectivity on patch.
tmp< scalarField > diffReflectivity(const label patchI, const label bandI=0, const vectorField *incomingDirection=nullptr, const scalarField *T=nullptr) const
Access boundary diffuse reflectivity on patch.
scalar faceEmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary emissivity on face.
A class for managing temporary objects.
Definition tmp.H:75
Namespace for radiation modelling.
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Vector< scalar > vector
Definition vector.H:57
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68