boundaryRadiationProperties.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-2020 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 \*---------------------------------------------------------------------------*/
27 
29 #include "radiationModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  namespace radiation
36  {
37  defineTypeNameAndDebug(boundaryRadiationProperties, 0);
38  }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const fvMesh& mesh
47 )
48 :
50  <
51  fvMesh,
54  >(mesh),
55  radBoundaryPropertiesPtrList_(mesh.boundary().size())
56 {
57  IOobject boundaryIO
58  (
59  boundaryRadiationProperties::typeName,
60  mesh.time().constant(),
61  mesh,
62  IOobject::MUST_READ,
63  IOobject::NO_WRITE,
64  false
65  );
66 
67  if (boundaryIO.typeHeaderOk<IOdictionary>(true))
68  {
69  const radiationModel& radiation =
70  mesh.lookupObject<radiationModel>
71  (
72  "radiationProperties"
73  );
74 
75  // Model number of bands
76  label nBands = radiation.nBands();
77 
78  const IOdictionary radiationDict(boundaryIO);
79 
80  forAll(mesh.boundary(), patchi)
81  {
82  const polyPatch& pp = mesh.boundaryMesh()[patchi];
83 
84  if (radiationDict.isDict(pp.name()))
85  {
86  const dictionary& dict = radiationDict.subDict(pp.name());
87 
88  radBoundaryPropertiesPtrList_[patchi].reset
89  (
91  );
92 
93  if (nBands != radBoundaryPropertiesPtrList_[patchi]->nBands())
94  {
96  << "Radiation bands : " << nBands << nl
97  << "Bands on patch : " << patchi << " is "
98  << radBoundaryPropertiesPtrList_[patchi]->nBands()
99  << abort(FatalError);
100  }
101  }
102  }
103  }
104 }
105 
106 
107 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
108 
111 (
112  const label patchi,
113  const label bandi,
114  vectorField* incomingDirection,
115  scalarField* T
116 ) const
117 {
118  if (radBoundaryPropertiesPtrList_[patchi])
119  {
120  return radBoundaryPropertiesPtrList_[patchi]->e
121  (
122  bandi,
123  incomingDirection,
124  T
125  );
126  }
127 
129  << "Patch : " << mesh().boundaryMesh()[patchi].name()
130  << " is not found in the boundaryRadiationProperties. "
131  << "Please add it"
132  << exit(FatalError);
133 
134  return tmp<scalarField>::New();
135 }
136 
137 
139 (
140  const label patchi,
141  const label faceI,
142  const label bandi,
143  vector incomingDirection,
144  scalar T
145 ) const
146 {
147  if (radBoundaryPropertiesPtrList_[patchi])
148  {
149  return radBoundaryPropertiesPtrList_[patchi]->e
150  (
151  faceI,
152  bandi,
153  incomingDirection,
154  T
155  );
156  }
157 
159  << "Patch : " << mesh().boundaryMesh()[patchi].name()
160  << " is not found in the boundaryRadiationProperties. "
161  << "Please add it"
162  << exit(FatalError);
163 
164  return Zero;
165 }
166 
167 
170 (
171  const label patchi,
172  const label bandi,
173  vectorField* incomingDirection,
174  scalarField* T
175 ) const
176 {
177  if (radBoundaryPropertiesPtrList_[patchi])
178  {
179  return radBoundaryPropertiesPtrList_[patchi]->a
180  (
181  bandi,
182  incomingDirection,
183  T
184  );
185  }
186 
188  << "Patch : " << mesh().boundaryMesh()[patchi].name()
189  << " is not found in the boundaryRadiationProperties. "
190  << "Please add it"
191  << exit(FatalError);
192 
193  return tmp<scalarField>::New();
194 }
195 
196 
198 (
199  const label patchi,
200  const label faceI,
201  const label bandi,
202  vector incomingDirection,
203  scalar T
204 ) const
205 {
206  if (radBoundaryPropertiesPtrList_[patchi])
207  {
208  return radBoundaryPropertiesPtrList_[patchi]->a
209  (
210  faceI,
211  bandi,
212  incomingDirection,
213  T
214  );
215  }
216 
218  << "Patch : " << mesh().boundaryMesh()[patchi].name()
219  << " is not found in the boundaryRadiationProperties. "
220  << "Please add it"
221  << exit(FatalError);
222 
223  return Zero;
224 }
225 
226 
229 (
230  const label patchi,
231  const label bandi,
232  vectorField* incomingDirection,
233  scalarField* T
234 ) const
235 {
236  if (radBoundaryPropertiesPtrList_[patchi])
237  {
238  return radBoundaryPropertiesPtrList_[patchi]->t
239  (
240  bandi,
241  incomingDirection,
242  T
243  );
244  }
245 
247  << "Patch : " << mesh().boundaryMesh()[patchi].name()
248  << " is not found in the boundaryRadiationProperties. "
249  << "Please add it"
250  << exit(FatalError);
251 
252  return tmp<scalarField>::New();
253 }
254 
255 
257 (
258  const label patchi,
259  const label faceI,
260  const label bandi,
261  vector incomingDirection,
262  scalar T
263 ) const
264 {
265  if (radBoundaryPropertiesPtrList_[patchi])
266  {
267  return radBoundaryPropertiesPtrList_[patchi]->t
268  (
269  faceI,
270  bandi,
271  incomingDirection,
272  T
273  );
274  }
275 
277  << "Patch : " << mesh().boundaryMesh()[patchi].name()
278  << " is not found in the boundaryRadiationProperties. "
279  << "Please add it"
280  << exit(FatalError);
281 
282  return Zero;
283 }
284 
285 
288 (
289  const label patchi,
290  const label bandi,
291  vectorField* incomingDirection,
292  scalarField* T
293 ) const
294 {
295  if (radBoundaryPropertiesPtrList_[patchi])
296  {
297  return radBoundaryPropertiesPtrList_[patchi]->rDiff
298  (
299  bandi,
300  incomingDirection,
301  T
302  );
303  }
304 
306  << "Patch : " << mesh().boundaryMesh()[patchi].name()
307  << " is not found in the boundaryRadiationProperties. "
308  << "Please add it"
309  << exit(FatalError);
310 
311  return tmp<scalarField>::New();
312 }
313 
314 
316 (
317  const label patchi,
318  const label faceI,
319  const label bandi,
320  vector incomingDirection,
321  scalar T
322 ) const
323 {
324  if (radBoundaryPropertiesPtrList_[patchi])
325  {
326  return radBoundaryPropertiesPtrList_[patchi]->rDiff
327  (
328  faceI,
329  bandi,
330  incomingDirection,
331  T
332  );
333  }
334 
336  << "Patch : " << mesh().boundaryMesh()[patchi].name()
337  << " is not found in the boundaryRadiationProperties. "
338  << "Please add it"
339  << exit(FatalError);
340 
341  return Zero;
342 }
343 
344 
347 (
348  const label patchi,
349  const label bandi,
350  vectorField* incomingDirection,
351  scalarField* T
352 ) const
353 {
354  if (radBoundaryPropertiesPtrList_[patchi])
355  {
356  return radBoundaryPropertiesPtrList_[patchi]->rSpec
357  (
358  bandi,
359  incomingDirection,
360  T
361  );
362  }
363 
365  << "Patch : " << mesh().boundaryMesh()[patchi].name()
366  << " is not found in the boundaryRadiationProperties. "
367  << "Please add it"
368  << exit(FatalError);
369 
370  return tmp<scalarField>::New();
371 }
372 
373 
375 (
376  const label patchi,
377  const label faceI,
378  const label bandi,
379  vector incomingDirection,
380  scalar T
381 ) const
382 {
383  if (radBoundaryPropertiesPtrList_[patchi])
384  {
385  return radBoundaryPropertiesPtrList_[patchi]->rSpec
386  (
387  faceI,
388  bandi,
389  incomingDirection,
390  T
391  );
392  }
393 
395  << "Patch : " << mesh().boundaryMesh()[patchi].name()
396  << " is not found in the boundaryRadiationProperties. "
397  << "Please add it"
398  << exit(FatalError);
399 
400  return Zero;
401 }
402 
403 
404 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::GeometricMeshObject
Definition: MeshObject.H:201
Foam::radiation::boundaryRadiationProperties::faceEmissivity
scalar faceEmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary emissivity on face.
Definition: boundaryRadiationProperties.C:139
Foam::radiation::boundaryRadiationProperties::boundaryRadiationProperties
boundaryRadiationProperties(const fvMesh &)
Construct given fvMesh.
Definition: boundaryRadiationProperties.C:45
Foam::radiation::boundaryRadiationProperties::faceSpecReflectivity
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.
Definition: boundaryRadiationProperties.C:375
Foam::radiation::boundaryRadiationProperties::specReflectivity
tmp< scalarField > specReflectivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary specular reflectivity on patch.
Definition: boundaryRadiationProperties.C:347
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::radiation::boundaryRadiationProperties::faceTransmissivity
scalar faceTransmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary transmissivity on face.
Definition: boundaryRadiationProperties.C:257
Foam::radiation::boundaryRadiationProperties::emissivity
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary emissivity on patch.
Definition: boundaryRadiationProperties.C:111
Foam::radiation::boundaryRadiationProperties::faceDiffReflectivity
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.
Definition: boundaryRadiationProperties.C:316
Foam::Field< vector >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Boundary radiation properties holder.
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::radiation::defineTypeNameAndDebug
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::dictionary::isDict
bool isDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Check if entry is found and is a sub-dictionary.
Definition: dictionaryI.H:147
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::radiation::boundaryRadiationProperties::faceAbsorptivity
scalar faceAbsorptivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary absorptivity on face.
Definition: boundaryRadiationProperties.C:198
Foam::Vector< scalar >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::radiation::boundaryRadiationProperties::absorptivity
tmp< scalarField > absorptivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary absorptivity on patch.
Definition: boundaryRadiationProperties.C:170
Foam::MeshObject
Templated abstract base-class for optional mesh objects used to automate their allocation to the mesh...
Definition: MeshObject.H:88
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
Foam::radiation::boundaryRadiationProperties::transmissivity
tmp< scalarField > transmissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary transmissivity on patch.
Definition: boundaryRadiationProperties.C:229
Foam::radiation::boundaryRadiationProperties::diffReflectivity
tmp< scalarField > diffReflectivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary diffuse reflectivity on patch.
Definition: boundaryRadiationProperties.C:288
radiationModel.H
boundaryRadiationProperties.H