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 -------------------------------------------------------------------------------
10 License
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 Class
27  Foam::radiation::boundaryRadiationProperties
28 
29 Description
30  Boundary radiation properties holder
31 
32 SourceFiles
33  boundaryRadiationProperties.C
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #ifndef boundaryRadiationProperties_H
38 #define boundaryRadiationProperties_H
39 
40 #include "MeshObject.H"
42 
43 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 
45 namespace Foam
46 {
47 
48 class fvMesh;
49 
50 namespace radiation
51 {
52 
53 /*---------------------------------------------------------------------------*\
54  Class boundaryRadiationProperties Declaration
55 \*---------------------------------------------------------------------------*/
56 
58 :
59  public MeshObject
60  <
61  fvMesh,
62  Foam::GeometricMeshObject,
63  boundaryRadiationProperties
64  >
65 {
66  // Private Data
67 
68  //- Per patch the boundaryRadiationProperties
70  radBoundaryPropertiesPtrList_;
71 
72  //- Per (face)zone the boundaryRadiationProperties
74  radZonePropertiesPtrList_;
75 
76 
77 public:
78 
79  // Declare name of the class and its debug switch
80  TypeName("boundaryRadiationProperties");
81 
82 
83  // Constructors
84 
85  //- Construct given fvMesh
86  explicit boundaryRadiationProperties(const fvMesh&);
87 
88 
89  // Member Functions
90 
91  //- Return identifiers of face zones activated for boundary radiation
92  const labelList faceZoneIds() const
93  {
94  DynamicList<label> ncZones;
95 
96  forAll(radZonePropertiesPtrList_, i)
97  {
98  if (radZonePropertiesPtrList_.test(i))
99  {
100  ncZones.append(i);
101  }
102  }
103 
104  return ncZones;
105  }
106 
107  //- Access boundary emissivity on patch
109  (
110  const label patchI,
111  const label bandI = 0,
112  vectorField* incomingDirection = nullptr,
113  scalarField* T = nullptr
114  ) const;
115 
116  //- Access boundary emissivity on face
117  scalar faceEmissivity
118  (
119  const label patchI,
120  const label faceI,
121  const label bandI = 0,
122  vector incomingDirection = Zero,
123  scalar T = 0
124  ) const;
125 
126  //- Access boundary absorptivity on patch
128  (
129  const label patchI,
130  const label bandI = 0,
131  vectorField* incomingDirection = nullptr,
132  scalarField* T = nullptr
133  ) const;
134 
135  //- Access boundary absorptivity on face
136  scalar faceAbsorptivity
137  (
138  const label patchI,
139  const label faceI,
140  const label bandI = 0,
141  vector incomingDirection = Zero,
142  scalar T = 0
143  ) const;
144 
145  //- Access boundary transmissivity on patch
147  (
148  const label patchI,
149  const label bandI = 0,
150  vectorField* incomingDirection = nullptr,
151  scalarField* T = nullptr
152  ) const;
153 
154  //- Access boundary transmissivity on face
155  scalar faceTransmissivity
156  (
157  const label patchI,
158  const label faceI,
159  const label bandI = 0,
160  vector incomingDirection = Zero,
161  scalar T = 0
162  ) const;
163 
164  // Specialisations for faceZones
165 
166  //- Access transmissivity on set of (internal) faces. Zone name only
167  // used to lookup the properties in boundaryRadiationProperties
169  (
170  const label zoneI,
171  const labelUList& faceIDs, // internal faces
172  const label bandI = 0,
173  vector incomingDirection = Zero,
174  scalar T = 0
175  ) const;
176 
177 
178  //- Access boundary diffuse reflectivity on patch
180  (
181  const label patchI,
182  const label bandI = 0,
183  vectorField* incomingDirection = nullptr,
184  scalarField* T = nullptr
185  ) const;
186 
187  //- Access boundary diffuse reflectivity on face
188  scalar faceDiffReflectivity
189  (
190  const label patchI,
191  const label faceI,
192  const label bandI = 0,
193  vector incomingDirection = Zero,
194  scalar T = 0
195  ) const;
196 
197  //- Access boundary specular reflectivity on patch
199  (
200  const label patchI,
201  const label bandI = 0,
202  vectorField* incomingDirection = nullptr,
203  scalarField* T = nullptr
204  ) const;
205 
206  //- Access boundary specular reflectivity on face
207  scalar faceSpecReflectivity
208  (
209  const label patchI,
210  const label faceI,
211  const label bandI = 0,
212  vector incomingDirection = Zero,
213  scalar T = 0
214  ) const;
215 
216 
217  //- Destructor
218  ~boundaryRadiationProperties() = default;
219 };
220 
221 
222 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
223 
224 } // End namespace radiation
225 } // End namespace Foam
226 
227 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
228 
229 
230 #endif
231 
232 // ************************************************************************* //
const T * test(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrListI.H:127
scalar faceAbsorptivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary absorptivity on face.
scalar faceTransmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary transmissivity on face.
scalar faceEmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary emissivity on face.
TypeName("boundaryRadiationProperties")
tmp< scalarField > specReflectivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary specular reflectivity on patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
tmp< scalarField > absorptivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary absorptivity on patch.
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:85
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
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.
boundaryRadiationProperties(const fvMesh &)
Construct given fvMesh.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
tmp< scalarField > diffReflectivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary diffuse reflectivity on patch.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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.
const labelList faceZoneIds() const
Return identifiers of face zones activated for boundary radiation.
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary transmissivity on patch.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
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.
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary emissivity on patch.
~boundaryRadiationProperties()=default
Destructor.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127