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,
65  false
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  const auto* ePtr = radiationDict.findEntry(pp.name());
93 
94  if (ePtr && ePtr->isDict())
95  {
96  radBoundaryPropertiesPtrList_.set
97  (
98  patchi,
100  );
101 
102  matchedEntries.insert(pp.name());
103 
104  if
105  (
106  nBands
107  != radBoundaryPropertiesPtrList_[patchi].nBands()
108  )
109  {
111  << "Radiation bands : " << nBands << nl
112  << "Bands on patch : " << patchi << " is "
113  << radBoundaryPropertiesPtrList_[patchi].nBands()
114  << abort(FatalError);
115  }
116  }
117  }
118  }
119 
120 
121  // Match faceZones if any dictionary entries have not been used for
122  // patch matching.
123  //
124  // Note: radiation properties are hardcoded to take patch reference.
125  // Supply patch0 for now.
126  {
127  const auto& dummyRef = mesh.boundaryMesh()[0];
128 
129  const auto& fzs = mesh.faceZones();
130 
131  for (const auto& fz : fzs)
132  {
133  const label zonei = fz.index();
134 
135  if (!matchedEntries.found(fz.name()))
136  {
137  // Note: avoid wildcard matches. Assume user explicitly
138  // provided information for faceZones.
139  const auto* ePtr = radiationDict.findEntry
140  (
141  fz.name(),
143  );
144 
145  // Skip face zones that are not used in boundary radiation
146  if (!ePtr)
147  {
148  continue;
149  }
150 
151  if (ePtr->isDict())
152  {
153  const dictionary& dict = ePtr->dict();
154 
155  radZonePropertiesPtrList_.set
156  (
157  zonei,
159  (
160  dict,
161  dummyRef
162  )
163  );
164 
165  matchedEntries.insert(fz.name());
166 
167  if
168  (
169  nBands
170  != radZonePropertiesPtrList_[zonei].nBands()
171  )
172  {
174  << "Radiation bands : " << nBands << nl
175  << "Bands on zone : " << zonei << " is "
176  << radBoundaryPropertiesPtrList_
177  [
178  zonei
179  ].nBands()
180  << abort(FatalError);
181  }
182  }
183  }
184  }
185  }
186  }
187 }
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
194 (
195  const label patchi,
196  const label bandi,
197  vectorField* incomingDirection,
198  scalarField* T
199 ) const
200 {
201  if (radBoundaryPropertiesPtrList_.set(patchi))
202  {
203  return radBoundaryPropertiesPtrList_[patchi].e
204  (
205  bandi,
206  incomingDirection,
207  T
208  );
209  }
210 
212  << "Patch : " << mesh().boundaryMesh()[patchi].name()
213  << " is not found in the boundaryRadiationProperties. "
214  << "Please add it"
216 
217  return tmp<scalarField>::New();
218 }
219 
220 
222 (
223  const label patchi,
224  const label facei,
225  const label bandi,
226  vector incomingDirection,
227  scalar T
228 ) const
229 {
230  if (radBoundaryPropertiesPtrList_.set(patchi))
231  {
232  return radBoundaryPropertiesPtrList_[patchi].e
233  (
234  facei,
235  bandi,
236  incomingDirection,
237  T
238  );
239  }
240 
242  << "Patch : " << mesh().boundaryMesh()[patchi].name()
243  << " is not found in the boundaryRadiationProperties. "
244  << "Please add it"
245  << exit(FatalError);
247  return Zero;
248 }
249 
250 
253 (
254  const label patchi,
255  const label bandi,
256  vectorField* incomingDirection,
257  scalarField* T
258 ) const
259 {
260  if (radBoundaryPropertiesPtrList_.set(patchi))
261  {
262  return radBoundaryPropertiesPtrList_[patchi].a
263  (
264  bandi,
265  incomingDirection,
266  T
267  );
268  }
269 
271  << "Patch : " << mesh().boundaryMesh()[patchi].name()
272  << " is not found in the boundaryRadiationProperties. "
273  << "Please add it"
275 
276  return tmp<scalarField>::New();
277 }
278 
279 
281 (
282  const label patchi,
283  const label facei,
284  const label bandi,
285  vector incomingDirection,
286  scalar T
287 ) const
288 {
289  if (radBoundaryPropertiesPtrList_.set(patchi))
290  {
291  return radBoundaryPropertiesPtrList_[patchi].a
292  (
293  facei,
294  bandi,
295  incomingDirection,
296  T
297  );
298  }
299 
301  << "Patch : " << mesh().boundaryMesh()[patchi].name()
302  << " is not found in the boundaryRadiationProperties. "
303  << "Please add it"
304  << exit(FatalError);
306  return Zero;
307 }
308 
309 
312 (
313  const label patchi,
314  const label bandi,
315  vectorField* incomingDirection,
316  scalarField* T
317 ) const
318 {
319  if (radBoundaryPropertiesPtrList_.set(patchi))
320  {
321  return radBoundaryPropertiesPtrList_[patchi].t
322  (
323  bandi,
324  incomingDirection,
325  T
326  );
327  }
328 
330  << "Patch : " << mesh().boundaryMesh()[patchi].name()
331  << " is not found in the boundaryRadiationProperties. "
332  << "Please add it"
334 
335  return tmp<scalarField>::New();
336 }
337 
338 
340 (
341  const label patchi,
342  const label facei,
343  const label bandi,
344  vector incomingDirection,
345  scalar T
346 ) const
347 {
348  if (radBoundaryPropertiesPtrList_.set(patchi))
349  {
350  return radBoundaryPropertiesPtrList_[patchi].t
351  (
352  facei,
353  bandi,
354  incomingDirection,
355  T
356  );
357  }
358 
360  << "Patch : " << mesh().boundaryMesh()[patchi].name()
361  << " is not found in the boundaryRadiationProperties. "
362  << "Please add it"
363  << exit(FatalError);
365  return Zero;
366 }
367 
368 
371 (
372  const label zonei,
373  const labelUList& faceIDs,
374  const label bandi,
375  vector incomingDirection,
376  scalar T
377 ) const
378 {
379  if (radZonePropertiesPtrList_.set(zonei))
380  {
381  auto tfld = tmp<scalarField>::New(faceIDs.size());
382  auto& fld = tfld.ref();
383  forAll(fld, i)
384  {
385  fld[i] = radZonePropertiesPtrList_[zonei].t
386  (
387  faceIDs[i],
388  bandi,
389  incomingDirection,
390  T
391  );
392  }
393  return tfld;
394  }
395 
397  << "Zone : " << mesh().faceZones()[zonei].name()
398  << " is not found in the boundaryRadiationProperties. "
399  << "Please add it"
400  << exit(FatalError);
402  return tmp<scalarField>::New();
403 }
404 
405 
408 (
409  const label patchi,
410  const label bandi,
411  vectorField* incomingDirection,
412  scalarField* T
413 ) const
414 {
415  if (radBoundaryPropertiesPtrList_.set(patchi))
416  {
417  return radBoundaryPropertiesPtrList_[patchi].rDiff
418  (
419  bandi,
420  incomingDirection,
421  T
422  );
423  }
424 
426  << "Patch : " << mesh().boundaryMesh()[patchi].name()
427  << " is not found in the boundaryRadiationProperties. "
428  << "Please add it"
430 
431  return tmp<scalarField>::New();
432 }
433 
434 
436 (
437  const label patchi,
438  const label facei,
439  const label bandi,
440  vector incomingDirection,
441  scalar T
442 ) const
443 {
444  if (radBoundaryPropertiesPtrList_.set(patchi))
445  {
446  return radBoundaryPropertiesPtrList_[patchi].rDiff
447  (
448  facei,
449  bandi,
450  incomingDirection,
451  T
452  );
453  }
454 
456  << "Patch : " << mesh().boundaryMesh()[patchi].name()
457  << " is not found in the boundaryRadiationProperties. "
458  << "Please add it"
459  << exit(FatalError);
461  return Zero;
462 }
463 
464 
467 (
468  const label patchi,
469  const label bandi,
470  vectorField* incomingDirection,
471  scalarField* T
472 ) const
473 {
474  if (radBoundaryPropertiesPtrList_.set(patchi))
475  {
476  return radBoundaryPropertiesPtrList_[patchi].rSpec
477  (
478  bandi,
479  incomingDirection,
480  T
481  );
482  }
483 
485  << "Patch : " << mesh().boundaryMesh()[patchi].name()
486  << " is not found in the boundaryRadiationProperties. "
487  << "Please add it"
489 
490  return tmp<scalarField>::New();
491 }
492 
493 
495 (
496  const label patchi,
497  const label facei,
498  const label bandi,
499  vector incomingDirection,
500  scalar T
501 ) const
502 {
503  if (radBoundaryPropertiesPtrList_.set(patchi))
504  {
505  return radBoundaryPropertiesPtrList_[patchi].rSpec
506  (
507  facei,
508  bandi,
509  incomingDirection,
510  T
511  );
512  }
513 
515  << "Patch : " << mesh().boundaryMesh()[patchi].name()
516  << " is not found in the boundaryRadiationProperties. "
517  << "Please add it"
518  << exit(FatalError);
519 
520  return Zero;
521 }
522 
523 
524 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:584
faceListList boundary
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
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:578
bool found(const Key &key) const
True if hashed key is found in table.
Definition: HashTableI.H:93
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
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:361
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:227
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:413
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:84
dynamicFvMesh & mesh
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:203
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:89
const faceZoneMesh & faceZones() const noexcept
Return face zone mesh.
Definition: polyMesh.H:646
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:79
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:80
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:166
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157