Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
26\*---------------------------------------------------------------------------*/
27
29#include "radiationModel.H"
31// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32
33namespace Foam
34{
35 namespace radiation
36 {
38 }
39}
40
41
42// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43
45(
46 const fvMesh& mesh
47)
48:
49 MeshObject_type(mesh),
50
51 radBoundaryPropertiesPtrList_(mesh.boundary().size()),
52 radZonePropertiesPtrList_(mesh.faceZones().size())
53{
54 IOobject boundaryIO
55 (
56 boundaryRadiationProperties::typeName,
57 mesh.time().constant(),
58 mesh,
62 );
63
64 if (boundaryIO.typeHeaderOk<IOdictionary>(true))
65 {
67 mesh.lookupObject<radiationModel>
68 (
69 "radiationProperties"
70 );
71
72 // Model number of bands
73 label nBands = radiation.nBands();
74
75 // Load in dictionary
76 const IOdictionary radiationDict(boundaryIO);
77
78
79 wordHashSet matchedEntries;
80
81 // Match patches
82 {
83 const auto& pbm = mesh.boundaryMesh();
84 for (const auto& pp : pbm)
85 {
86 const label patchi = pp.index();
87
88 // Allow wildcard matches
89 const auto* eptr =
90 radiationDict.findEntry(pp.name(), keyType::REGEX);
91
92 if (eptr && eptr->isDict())
93 {
94 radBoundaryPropertiesPtrList_.set
95 (
96 patchi,
98 );
99
100 matchedEntries.insert(pp.name());
101
102 if
103 (
104 nBands
105 != radBoundaryPropertiesPtrList_[patchi].nBands()
106 )
107 {
109 << "Radiation bands : " << nBands << nl
110 << "Bands on patch : " << patchi << " is "
111 << radBoundaryPropertiesPtrList_[patchi].nBands()
112 << abort(FatalError);
113 }
114 }
115 }
116 }
117
118
119 // Match faceZones if any dictionary entries have not been used for
120 // patch matching.
121 //
122 // Note: radiation properties are hardcoded to take patch reference.
123 // Supply patch0 for now.
124 {
125 const auto& dummyRef = mesh.boundaryMesh()[0];
126
127 const auto& fzs = mesh.faceZones();
128
129 for (const auto& fz : fzs)
130 {
131 const label zonei = fz.index();
132
133 if (!matchedEntries.contains(fz.name()))
134 {
135 // Note: avoid wildcard matches. Assume user explicitly
136 // provided information for faceZones.
137
138 const auto* eptr = radiationDict.findEntry
139 (
140 fz.name(),
142 );
143
144 if (eptr && eptr->isDict())
145 {
146 radZonePropertiesPtrList_.set
147 (
148 zonei,
150 (
151 eptr->dict(),
152 dummyRef
153 )
154 );
155
156 matchedEntries.insert(fz.name());
157
158 if
159 (
160 nBands
161 != radZonePropertiesPtrList_[zonei].nBands()
162 )
163 {
165 << "Radiation bands : " << nBands << nl
166 << "Bands on zone : " << zonei << " is "
167 << radBoundaryPropertiesPtrList_
168 [
169 zonei
170 ].nBands()
171 << abort(FatalError);
172 }
173 }
174 }
175 }
176 }
178}
179
180
181// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
182
185(
186 const label patchi,
187 const label bandi,
188 const vectorField* incomingDirection,
189 const scalarField* T
190) const
191{
192 if (radBoundaryPropertiesPtrList_.set(patchi))
193 {
194 return radBoundaryPropertiesPtrList_[patchi].e
195 (
196 bandi,
197 incomingDirection,
198 T
199 );
200 }
201
203 << "Patch : " << mesh().boundaryMesh()[patchi].name()
204 << " is not found in the boundaryRadiationProperties. "
205 << "Please add it"
206 << exit(FatalError);
207
208 return nullptr;
209}
210
211
213(
214 const label patchi,
215 const label facei,
216 const label bandi,
217 vector incomingDirection,
218 scalar T
219) const
220{
221 if (radBoundaryPropertiesPtrList_.set(patchi))
222 {
223 return radBoundaryPropertiesPtrList_[patchi].e
224 (
225 facei,
226 bandi,
227 incomingDirection,
228 T
229 );
230 }
231
233 << "Patch : " << mesh().boundaryMesh()[patchi].name()
234 << " is not found in the boundaryRadiationProperties. "
235 << "Please add it"
237
238 return Zero;
239}
240
241
244(
245 const label patchi,
246 const label bandi,
247 const vectorField* incomingDirection,
248 const scalarField* T
249) const
250{
251 if (radBoundaryPropertiesPtrList_.set(patchi))
252 {
253 return radBoundaryPropertiesPtrList_[patchi].a
254 (
255 bandi,
256 incomingDirection,
257 T
258 );
259 }
260
262 << "Patch : " << mesh().boundaryMesh()[patchi].name()
263 << " is not found in the boundaryRadiationProperties. "
264 << "Please add it"
265 << exit(FatalError);
266
267 return nullptr;
268}
269
270
272(
273 const label patchi,
274 const label facei,
275 const label bandi,
276 vector incomingDirection,
277 scalar T
278) const
279{
280 if (radBoundaryPropertiesPtrList_.set(patchi))
281 {
282 return radBoundaryPropertiesPtrList_[patchi].a
283 (
284 facei,
285 bandi,
286 incomingDirection,
287 T
288 );
289 }
290
292 << "Patch : " << mesh().boundaryMesh()[patchi].name()
293 << " is not found in the boundaryRadiationProperties. "
294 << "Please add it"
296
297 return Zero;
298}
299
300
303(
304 const label patchi,
305 const label bandi,
306 const vectorField* incomingDirection,
307 const scalarField* T
308) const
309{
310 if (radBoundaryPropertiesPtrList_.set(patchi))
311 {
312 return radBoundaryPropertiesPtrList_[patchi].t
313 (
314 bandi,
315 incomingDirection,
316 T
317 );
318 }
319
321 << "Patch : " << mesh().boundaryMesh()[patchi].name()
322 << " is not found in the boundaryRadiationProperties. "
323 << "Please add it"
324 << exit(FatalError);
325
326 return nullptr;
327}
328
329
331(
332 const label patchi,
333 const label facei,
334 const label bandi,
335 vector incomingDirection,
336 scalar T
337) const
338{
339 if (radBoundaryPropertiesPtrList_.set(patchi))
340 {
341 return radBoundaryPropertiesPtrList_[patchi].t
342 (
343 facei,
344 bandi,
345 incomingDirection,
346 T
347 );
348 }
349
351 << "Patch : " << mesh().boundaryMesh()[patchi].name()
352 << " is not found in the boundaryRadiationProperties. "
353 << "Please add it"
355
356 return Zero;
357}
358
359
362(
363 const label zonei,
364 const labelUList& faceIDs,
365 const label bandi,
366 vector incomingDirection,
367 scalar T
368) const
369{
370 if (radZonePropertiesPtrList_.set(zonei))
371 {
372 auto tfld = tmp<scalarField>::New(faceIDs.size());
373 auto& fld = tfld.ref();
374 forAll(fld, i)
375 {
376 fld[i] = radZonePropertiesPtrList_[zonei].t
377 (
378 faceIDs[i],
379 bandi,
380 incomingDirection,
381 T
382 );
383 }
384 return tfld;
385 }
386
388 << "Zone : " << mesh().faceZones()[zonei].name()
389 << " is not found in the boundaryRadiationProperties. "
390 << "Please add it"
392
393 return nullptr;
394}
395
396
399(
400 const label patchi,
401 const label bandi,
402 const vectorField* incomingDirection,
403 const scalarField* T
404) const
405{
406 if (radBoundaryPropertiesPtrList_.set(patchi))
407 {
408 return radBoundaryPropertiesPtrList_[patchi].rDiff
409 (
410 bandi,
411 incomingDirection,
412 T
413 );
414 }
415
417 << "Patch : " << mesh().boundaryMesh()[patchi].name()
418 << " is not found in the boundaryRadiationProperties. "
419 << "Please add it"
420 << exit(FatalError);
421
422 return nullptr;
423}
424
425
427(
428 const label patchi,
429 const label facei,
430 const label bandi,
431 vector incomingDirection,
432 scalar T
433) const
434{
435 if (radBoundaryPropertiesPtrList_.set(patchi))
436 {
437 return radBoundaryPropertiesPtrList_[patchi].rDiff
438 (
439 facei,
440 bandi,
441 incomingDirection,
442 T
443 );
444 }
445
447 << "Patch : " << mesh().boundaryMesh()[patchi].name()
448 << " is not found in the boundaryRadiationProperties. "
449 << "Please add it"
451
452 return Zero;
453}
454
455
458(
459 const label patchi,
460 const label bandi,
461 const vectorField* incomingDirection,
462 const scalarField* T
463) const
464{
465 if (radBoundaryPropertiesPtrList_.set(patchi))
466 {
467 return radBoundaryPropertiesPtrList_[patchi].rSpec
468 (
469 bandi,
470 incomingDirection,
471 T
472 );
473 }
474
476 << "Patch : " << mesh().boundaryMesh()[patchi].name()
477 << " is not found in the boundaryRadiationProperties. "
478 << "Please add it"
479 << exit(FatalError);
480
481 return nullptr;
482}
483
484
486(
487 const label patchi,
488 const label facei,
489 const label bandi,
490 vector incomingDirection,
491 scalar T
492) const
493{
494 if (radBoundaryPropertiesPtrList_.set(patchi))
495 {
496 return radBoundaryPropertiesPtrList_[patchi].rSpec
497 (
498 facei,
499 bandi,
500 incomingDirection,
501 T
502 );
503 }
504
506 << "Patch : " << mesh().boundaryMesh()[patchi].name()
507 << " is not found in the boundaryRadiationProperties. "
508 << "Please add it"
509 << exit(FatalError);
510
511 return Zero;
512}
513
514
515// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
bool contains(const Key &key) const
True if hashed key is contained (found) in table.
Definition HashTableI.H:72
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info. A void type suppresses trait and t...
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
@ LITERAL
String literal.
Definition keyType.H:82
@ REGEX
Regular expression.
Definition keyType.H:83
static autoPtr< boundaryRadiationPropertiesPatch > New(const dictionary &dict, const polyPatch &pp)
Selector.
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 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.
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.
Top level model for radiation modelling.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
faceListList boundary
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
Namespace for radiation modelling.
Namespace for OpenFOAM.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299