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-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 \*---------------------------------------------------------------------------*/
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  radZonePropertiesPtrList_(mesh.faceZones().size())
57 {
58  IOobject boundaryIO
59  (
60  boundaryRadiationProperties::typeName,
61  mesh.time().constant(),
62  mesh,
66  );
67 
68  if (boundaryIO.typeHeaderOk<IOdictionary>(true))
69  {
70  const radiationModel& radiation =
72  (
73  "radiationProperties"
74  );
75 
76  // Model number of bands
77  label nBands = radiation.nBands();
78 
79  // Load in dictionary
80  const IOdictionary radiationDict(boundaryIO);
81 
82 
83  wordHashSet matchedEntries;
84 
85  // Match patches
86  {
87  const auto& pbm = mesh.boundaryMesh();
88  for (const auto& pp : pbm)
89  {
90  const label patchi = pp.index();
91 
92  // Allow wildcard matches
93  const auto* eptr =
94  radiationDict.findEntry(pp.name(), keyType::REGEX);
95 
96  if (eptr && eptr->isDict())
97  {
98  radBoundaryPropertiesPtrList_.set
99  (
100  patchi,
102  );
103 
104  matchedEntries.insert(pp.name());
105 
106  if
107  (
108  nBands
109  != radBoundaryPropertiesPtrList_[patchi].nBands()
110  )
111  {
113  << "Radiation bands : " << nBands << nl
114  << "Bands on patch : " << patchi << " is "
115  << radBoundaryPropertiesPtrList_[patchi].nBands()
116  << abort(FatalError);
117  }
118  }
119  }
120  }
121 
122 
123  // Match faceZones if any dictionary entries have not been used for
124  // patch matching.
125  //
126  // Note: radiation properties are hardcoded to take patch reference.
127  // Supply patch0 for now.
128  {
129  const auto& dummyRef = mesh.boundaryMesh()[0];
130 
131  const auto& fzs = mesh.faceZones();
132 
133  for (const auto& fz : fzs)
134  {
135  const label zonei = fz.index();
136 
137  if (!matchedEntries.contains(fz.name()))
138  {
139  // Note: avoid wildcard matches. Assume user explicitly
140  // provided information for faceZones.
141 
142  const auto* eptr = radiationDict.findEntry
143  (
144  fz.name(),
146  );
147 
148  if (eptr && eptr->isDict())
149  {
150  radZonePropertiesPtrList_.set
151  (
152  zonei,
154  (
155  eptr->dict(),
156  dummyRef
157  )
158  );
159 
160  matchedEntries.insert(fz.name());
161 
162  if
163  (
164  nBands
165  != radZonePropertiesPtrList_[zonei].nBands()
166  )
167  {
169  << "Radiation bands : " << nBands << nl
170  << "Bands on zone : " << zonei << " is "
171  << radBoundaryPropertiesPtrList_
172  [
173  zonei
174  ].nBands()
175  << abort(FatalError);
176  }
177  }
178  }
179  }
180  }
181  }
182 }
183 
184 
185 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
186 
189 (
190  const label patchi,
191  const label bandi,
192  vectorField* incomingDirection,
193  scalarField* T
194 ) const
195 {
196  if (radBoundaryPropertiesPtrList_.set(patchi))
197  {
198  return radBoundaryPropertiesPtrList_[patchi].e
199  (
200  bandi,
201  incomingDirection,
202  T
203  );
204  }
205 
207  << "Patch : " << mesh().boundaryMesh()[patchi].name()
208  << " is not found in the boundaryRadiationProperties. "
209  << "Please add it"
211 
212  return tmp<scalarField>::New();
213 }
214 
215 
217 (
218  const label patchi,
219  const label facei,
220  const label bandi,
221  vector incomingDirection,
222  scalar T
223 ) const
224 {
225  if (radBoundaryPropertiesPtrList_.set(patchi))
226  {
227  return radBoundaryPropertiesPtrList_[patchi].e
228  (
229  facei,
230  bandi,
231  incomingDirection,
232  T
233  );
234  }
235 
237  << "Patch : " << mesh().boundaryMesh()[patchi].name()
238  << " is not found in the boundaryRadiationProperties. "
239  << "Please add it"
240  << exit(FatalError);
242  return Zero;
243 }
244 
245 
248 (
249  const label patchi,
250  const label bandi,
251  vectorField* incomingDirection,
252  scalarField* T
253 ) const
254 {
255  if (radBoundaryPropertiesPtrList_.set(patchi))
256  {
257  return radBoundaryPropertiesPtrList_[patchi].a
258  (
259  bandi,
260  incomingDirection,
261  T
262  );
263  }
264 
266  << "Patch : " << mesh().boundaryMesh()[patchi].name()
267  << " is not found in the boundaryRadiationProperties. "
268  << "Please add it"
270 
271  return tmp<scalarField>::New();
272 }
273 
274 
276 (
277  const label patchi,
278  const label facei,
279  const label bandi,
280  vector incomingDirection,
281  scalar T
282 ) const
283 {
284  if (radBoundaryPropertiesPtrList_.set(patchi))
285  {
286  return radBoundaryPropertiesPtrList_[patchi].a
287  (
288  facei,
289  bandi,
290  incomingDirection,
291  T
292  );
293  }
294 
296  << "Patch : " << mesh().boundaryMesh()[patchi].name()
297  << " is not found in the boundaryRadiationProperties. "
298  << "Please add it"
299  << exit(FatalError);
301  return Zero;
302 }
303 
304 
307 (
308  const label patchi,
309  const label bandi,
310  vectorField* incomingDirection,
311  scalarField* T
312 ) const
313 {
314  if (radBoundaryPropertiesPtrList_.set(patchi))
315  {
316  return radBoundaryPropertiesPtrList_[patchi].t
317  (
318  bandi,
319  incomingDirection,
320  T
321  );
322  }
323 
325  << "Patch : " << mesh().boundaryMesh()[patchi].name()
326  << " is not found in the boundaryRadiationProperties. "
327  << "Please add it"
329 
330  return tmp<scalarField>::New();
331 }
332 
333 
335 (
336  const label patchi,
337  const label facei,
338  const label bandi,
339  vector incomingDirection,
340  scalar T
341 ) const
342 {
343  if (radBoundaryPropertiesPtrList_.set(patchi))
344  {
345  return radBoundaryPropertiesPtrList_[patchi].t
346  (
347  facei,
348  bandi,
349  incomingDirection,
350  T
351  );
352  }
353 
355  << "Patch : " << mesh().boundaryMesh()[patchi].name()
356  << " is not found in the boundaryRadiationProperties. "
357  << "Please add it"
358  << exit(FatalError);
360  return Zero;
361 }
362 
363 
366 (
367  const label zonei,
368  const labelUList& faceIDs,
369  const label bandi,
370  vector incomingDirection,
371  scalar T
372 ) const
373 {
374  if (radZonePropertiesPtrList_.set(zonei))
375  {
376  auto tfld = tmp<scalarField>::New(faceIDs.size());
377  auto& fld = tfld.ref();
378  forAll(fld, i)
379  {
380  fld[i] = radZonePropertiesPtrList_[zonei].t
381  (
382  faceIDs[i],
383  bandi,
384  incomingDirection,
385  T
386  );
387  }
388  return tfld;
389  }
390 
392  << "Zone : " << mesh().faceZones()[zonei].name()
393  << " is not found in the boundaryRadiationProperties. "
394  << "Please add it"
395  << exit(FatalError);
397  return tmp<scalarField>::New();
398 }
399 
400 
403 (
404  const label patchi,
405  const label bandi,
406  vectorField* incomingDirection,
407  scalarField* T
408 ) const
409 {
410  if (radBoundaryPropertiesPtrList_.set(patchi))
411  {
412  return radBoundaryPropertiesPtrList_[patchi].rDiff
413  (
414  bandi,
415  incomingDirection,
416  T
417  );
418  }
419 
421  << "Patch : " << mesh().boundaryMesh()[patchi].name()
422  << " is not found in the boundaryRadiationProperties. "
423  << "Please add it"
425 
426  return tmp<scalarField>::New();
427 }
428 
429 
431 (
432  const label patchi,
433  const label facei,
434  const label bandi,
435  vector incomingDirection,
436  scalar T
437 ) const
438 {
439  if (radBoundaryPropertiesPtrList_.set(patchi))
440  {
441  return radBoundaryPropertiesPtrList_[patchi].rDiff
442  (
443  facei,
444  bandi,
445  incomingDirection,
446  T
447  );
448  }
449 
451  << "Patch : " << mesh().boundaryMesh()[patchi].name()
452  << " is not found in the boundaryRadiationProperties. "
453  << "Please add it"
454  << exit(FatalError);
456  return Zero;
457 }
458 
459 
462 (
463  const label patchi,
464  const label bandi,
465  vectorField* incomingDirection,
466  scalarField* T
467 ) const
468 {
469  if (radBoundaryPropertiesPtrList_.set(patchi))
470  {
471  return radBoundaryPropertiesPtrList_[patchi].rSpec
472  (
473  bandi,
474  incomingDirection,
475  T
476  );
477  }
478 
480  << "Patch : " << mesh().boundaryMesh()[patchi].name()
481  << " is not found in the boundaryRadiationProperties. "
482  << "Please add it"
484 
485  return tmp<scalarField>::New();
486 }
487 
488 
490 (
491  const label patchi,
492  const label facei,
493  const label bandi,
494  vector incomingDirection,
495  scalar T
496 ) const
497 {
498  if (radBoundaryPropertiesPtrList_.set(patchi))
499  {
500  return radBoundaryPropertiesPtrList_[patchi].rSpec
501  (
502  facei,
503  bandi,
504  incomingDirection,
505  T
506  );
507  }
508 
510  << "Patch : " << mesh().boundaryMesh()[patchi].name()
511  << " is not found in the boundaryRadiationProperties. "
512  << "Please add it"
513  << exit(FatalError);
514 
515  return Zero;
516 }
517 
518 
519 // ************************************************************************* //
faceListList boundary
const polyBoundaryMesh & pbm
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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...
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
scalar faceEmissivity(const label patchI, const label faceI, const label bandI=0, vector incomingDirection=Zero, scalar T=0) const
Access boundary emissivity on face.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
defineTypeNameAndDebug(cloudAbsorptionEmission, 0)
tmp< scalarField > specReflectivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary specular reflectivity on patch.
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:232
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
#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
dynamicFvMesh & mesh
bool contains(const Key &key) const
True if hashed key is contained (found) in table.
Definition: HashTableI.H:72
const fvMesh & mesh() const noexcept
Reference to the mesh.
Definition: MeshObject.H:157
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
Top level model for radiation modelling.
String literal.
Definition: keyType.H:82
errorManip< error > abort(error &err)
Definition: errorManip.H:139
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...
const word & constant() const noexcept
Return constant name.
Definition: TimePathsI.H:112
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:670
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)
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
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 > transmissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary transmissivity on patch.
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.
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:84
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary emissivity on patch.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static autoPtr< boundaryRadiationPropertiesPatch > New(const dictionary &dict, const polyPatch &pp)
Selector.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Regular expression.
Definition: keyType.H:83
Do not request registration (bool: false)
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127