Loading...
Searching...
No Matches
PDRobstacle.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) 2016 Shell Research Ltd.
9 Copyright (C) 2019-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::PDRobstacle
29
30Description
31 Obstacle definitions for PDR
32
33SourceFiles
34 PDRobstacle.C
35 PDRobstacleIO.C
36 PDRobstacleRead.C
37
38\*---------------------------------------------------------------------------*/
39
40#ifndef Foam_PDRobstacle_H
41#define Foam_PDRobstacle_H
42
43#include "InfoProxy.H"
44#include "labelPair.H"
45#include "MeshedSurface.H"
46#include "MeshedSurfacesFwd.H"
47#include "boundBox.H"
48#include "DynamicList.H"
49#include "pointField.H"
50#include "volumeType.H"
52
53// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54
55namespace Foam
56{
58// Forward Declarations
59class boundBox;
60class PDRobstacle;
61
64
65namespace vtk
66{
67 class surfaceWriter;
68}
69
71/*---------------------------------------------------------------------------*\
72 Class PDRobstacle Declaration
73\*---------------------------------------------------------------------------*/
74
75class PDRobstacle
76{
77public:
78
79 //- Obstacle types (legacy numbering)
81 {
82 NONE = 0,
87 CUBOID = 6,
95 IGNITION = 41,
96 MESH_PLANE = 46,
97 IGNORE = 200
98 };
99
100
101 // Static Data Members
102
103 //- The max blowoff pressure [bar]
104 // Primarily to catch accidental input in Pa or mbar
105 static constexpr int maxBlowoffPressure = 10;
106
107
108 // Data Members
109
110 //- The group-id
111 label groupId;
112
113 //- The obstacle type-id
114 int typeId;
115
116 //- The x/y/z orientation (0,1,2)
119 //- Bias for position sorting
120 scalar sortBias;
121
122 //- The obstacle location.
123 // Lower corner for boxes, end-centre for cylinders
124 point pt;
125
126 //- The obstacle dimensions (for boxes)
127 vector span;
129 // Accessors for cylinders and diagonal blocks
130
131 inline scalar dia() const { return span[vector::X]; }
132 inline scalar theta() const { return span[vector::Y]; }
133 inline scalar len() const { return span[vector::Z]; }
134
135 inline scalar& dia() { return span[vector::X]; }
136 inline scalar& theta() { return span[vector::Y]; }
137 inline scalar& len() { return span[vector::Z]; }
138
139 union
141 scalar wa;
142 scalar slat_width;
143 scalar blowoff_press;
144 };
145 union
147 scalar wb;
149 };
150 scalar vbkge;
151 scalar xbkge;
152 scalar ybkge;
153 scalar zbkge;
155 union
157 int blowoff_type;
158 int inlet_dirn;
159 };
162
163public:
165 // Constructors
167 //- Construct zero-initialized
168 PDRobstacle();
169
170 //- Read construct as named dictionary
171 explicit PDRobstacle(Istream& is);
172
173
174 // Member Function Selectors
175
177 (
178 void,
180 read,
182 (
184 const dictionary& dict
185 ),
186 (obs, dict)
187 );
189
190 // Static Member Functions
191
192 //- Read obstacle files and add to the lists
193 // \return the total volume
194 static scalar legacyReadFiles
195 (
196 const fileName& obsFileDir,
197 const wordList& obsFileNames,
198 const boundBox& meshBb,
199 // output
201 DynamicList<PDRobstacle>& cylinders
202 );
203
204 //- Read obstacle files and set the lists
205 // \return the total volume
206 static scalar readFiles
207 (
208 const fileName& obsFileDir,
209 const wordList& obsFileNames,
210 const boundBox& meshBb,
211 // output
213 DynamicList<PDRobstacle>& cylinders
214 );
215
216
217 // Member Functions
218
219 //- Read name / dictionary
220 bool read(Istream& is);
221
222 //- Read the 'name' identifier if present
223 void readProperties(const dictionary& dict);
224
225 //- Obstacle position accessors
226 inline scalar x() const { return pt.x(); }
227 inline scalar y() const { return pt.y(); }
228 inline scalar z() const { return pt.z(); }
229 inline scalar& x() { return pt.x(); }
230 inline scalar& y() { return pt.y(); }
231 inline scalar& z() { return pt.z(); }
232
233
234 //- Is obstacle type id cylinder-like?
235 inline static bool isCylinder(const label id);
236
237 //- Is obstacle cylinder-like?
238 inline bool isCylinder() const;
239
240 //- Reset to a zero obstacle
241 void clear();
242
243 //- Scale obstacle dimensions by specified scaling factor
244 // Zero and negative factors are ignored
245 void scale(const scalar factor);
246
247 //- Volume of the obstacle
248 scalar volume() const;
249
250 //- True if the obstacle is considered to be too small
251 bool tooSmall(const scalar minWidth) const;
252
253 //- Set values from single-line, multi-column format.
254 // The only input format, but termed 'legacy' since it may
255 // be replaced in the near future.
256 // \return false if the scanning failed or if the obstacle type
257 // is not supported (or no longer supported)
260 const int groupTypeId,
261 const string& buffer,
262 const label lineNo = -1,
263 const word& inputFile = word::null
264 );
265
266 //- Trim obstacle to ensure it is within the specified bounding box
267 //- and return the intersection type.
268 // Returns UNKNOWN for unknown types and invalid bounding boxes
269 volumeType trim(const boundBox& bb);
270
271 //- Surface (points, faces) representation
272 meshedSurface surface() const;
273
274 //- Add pieces to vtp output
275 static label addPieces
277 vtk::surfaceWriter& surfWriter,
278 const UList<PDRobstacle>& list,
279 label pieceId = 0
280 );
281
282 //- Generate multi-piece VTK (vtp) file of obstacles
283 static void generateVtk
284 (
285 const fileName& outputDir,
286 const UList<PDRobstacle>& obslist,
287 const UList<PDRobstacle>& cyllist
288 );
289
290
291 // IOstream Operators
292
293 //- Return info proxy,
294 //- used to print information to a stream
296 {
297 return *this;
298 }
299
300 friend Istream& operator>>(Istream& is, PDRobstacle& obs);
301};
302
304// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305
306// Global Operators
307
308//- Compare according to x0 position
309bool operator<(const PDRobstacle& a, const PDRobstacle& b);
310
311//- For list output, assert that no obstacles are identical
312inline bool operator!=(const PDRobstacle& a, const PDRobstacle& b)
313{
314 return true;
315}
316
318// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319
320namespace PDRlegacy
321{
323/*---------------------------------------------------------------------------*\
324 Class obstacleGrouping Declaration
325\*---------------------------------------------------------------------------*/
326
327//- Locations for each instance of an obstacle group.
328class obstacleGrouping
329:
330 public DynamicList<point>
331{
332 //- Number of obstacles counted
333 label nObstacle_;
334
335 //- Number of cylinder-like obstacles counted
336 label nCylinder_;
338
339public:
340
341 //- Construct null
343 :
344 nObstacle_(0),
345 nCylinder_(0)
346 {}
347
348 //- Construct with one location (instance)
349 explicit obstacleGrouping(const vector& origin)
350 :
352 {
353 append(origin);
354 }
355
356 //- Clear obstacle count and locations
357 void clear()
358 {
359 nObstacle_ = 0;
360 nCylinder_ = 0;
362 }
363
364 //- Increment the number of obstacles
365 void addObstacle()
366 {
367 ++nObstacle_;
368 }
369
370 //- Increment the number of cylinder-like obstacles
371 void addCylinder()
373 ++nCylinder_;
374 }
375
376 //- The number of obstacles
377 label nObstacle() const
378 {
379 return nObstacle_;
381
382 //- The number of cylinder-like obstacles
383 label nCylinder() const
384 {
385 return nCylinder_;
386 }
387
388 //- The number of locations x number of obstacles
389 label nTotalObstacle() const
391 return size() * nObstacle_;
392 }
393
394 //- The number of locations x number of cylinder-like obstacles
395 label nTotalCylinder() const
396 {
397 return size() * nCylinder_;
398 }
399
400 //- The number of locations x number of obstacles
401 label nTotal() const
402 {
403 return size() * (nObstacle_ + nCylinder_);
404 }
405
406 //- Add a location
408
409 //- Add a location
410 void append(const scalar x, const scalar y, const scalar z)
411 {
412 append(point(x, y, z));
413 }
414};
415
416
417// Service Functions
418
419//- Read obstacle files, do counting only.
420// \return nObstacle, nCylinder read
422(
423 const fileName& obsFileDir,
424 const wordList& obsFileNames,
426);
427
428
429//- Read obstacle files and add to the lists
430// \return the total volume
432(
433 const fileName& obsFileDir,
434 const wordList& obsFileNames,
435 const Map<obstacleGrouping>& groups,
436 const boundBox& meshBb,
437 // output
440);
441
442
443// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
444
445} // End namespace PDRlegacy
446
448// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
449
450// Global Operators
451
452//- Locations for each instance of an obstacle group.
453inline Ostream& operator<<
454(
456 const PDRlegacy::obstacleGrouping& group
457)
458{
460 << group.size() << token::SPACE
461 << group.nObstacle() << token::SPACE
462 << group.nCylinder() << token::END_LIST;
464 return os;
465}
466
467
468// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
469
470} // End namespace Foam
472// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
473
474#include "PDRobstacleI.H"
475
476// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477
478#endif
480// ************************************************************************* //
scalar y
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
constexpr DynamicList() noexcept
A helper class for outputting values to Ostream.
Definition InfoProxy.H:49
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition Istream.H:60
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Locations for each instance of an obstacle group.
void addCylinder()
Increment the number of cylinder-like obstacles.
label nObstacle() const
The number of obstacles.
void addObstacle()
Increment the number of obstacles.
label nTotalObstacle() const
The number of locations x number of obstacles.
label nCylinder() const
The number of cylinder-like obstacles.
void append(const scalar x, const scalar y, const scalar z)
Add a location.
obstacleGrouping(const vector &origin)
Construct with one location (instance).
label nTotal() const
The number of locations x number of obstacles.
void clear()
Clear obstacle count and locations.
label nTotalCylinder() const
The number of locations x number of cylinder-like obstacles.
Obstacle definitions for PDR.
Definition PDRobstacle.H:71
scalar volume() const
Volume of the obstacle.
static label addPieces(vtk::surfaceWriter &surfWriter, const UList< PDRobstacle > &list, label pieceId=0)
Add pieces to vtp output.
scalar sortBias
Bias for position sorting.
bool isCylinder() const
Is obstacle cylinder-like?
point pt
The obstacle location.
PDRobstacle(Istream &is)
Read construct as named dictionary.
direction orient
The x/y/z orientation (0,1,2).
scalar x() const
Obstacle position accessors.
scalar y() const
scalar z() const
bool setFromLegacy(const int groupTypeId, const string &buffer, const label lineNo=-1, const word &inputFile=word::null)
Set values from single-line, multi-column format.
volumeType trim(const boundBox &bb)
Trim obstacle to ensure it is within the specified bounding box and return the intersection type.
scalar theta() const
scalar len() const
declareMemberFunctionSelectionTable(void, PDRobstacle, read, dictionary,(PDRobstacle &obs, const dictionary &dict),(obs, dict))
static void generateVtk(const fileName &outputDir, const UList< PDRobstacle > &obslist, const UList< PDRobstacle > &cyllist)
Generate multi-piece VTK (vtp) file of obstacles.
static scalar readFiles(const fileName &obsFileDir, const wordList &obsFileNames, const boundBox &meshBb, DynamicList< PDRobstacle > &blocks, DynamicList< PDRobstacle > &cylinders)
Read obstacle files and set the lists.
bool tooSmall(const scalar minWidth) const
True if the obstacle is considered to be too small.
meshedSurface surface() const
Surface (points, faces) representation.
label groupId
The group-id.
void readProperties(const dictionary &dict)
Read the 'name' identifier if present.
friend Istream & operator>>(Istream &is, PDRobstacle &obs)
void clear()
Reset to a zero obstacle.
scalar dia() const
static scalar legacyReadFiles(const fileName &obsFileDir, const wordList &obsFileNames, const boundBox &meshBb, DynamicList< PDRobstacle > &blocks, DynamicList< PDRobstacle > &cylinders)
Read obstacle files and add to the lists.
static constexpr int maxBlowoffPressure
The max blowoff pressure [bar].
vector span
The obstacle dimensions (for boxes).
bool read(Istream &is)
Read name / dictionary.
legacyTypes
Obstacle types (legacy numbering).
Definition PDRobstacle.H:78
@ OLD_BLOWOFF
ignored (old)
Definition PDRobstacle.H:88
@ IGNITION
ignored (old)
Definition PDRobstacle.H:92
@ OLD_INLET
ignored (old)
Definition PDRobstacle.H:87
@ NONE
Placeholder.
Definition PDRobstacle.H:79
int typeId
The obstacle type-id.
void scale(const scalar factor)
Scale obstacle dimensions by specified scaling factor.
InfoProxy< PDRobstacle > info() const noexcept
Return info proxy, used to print information to a stream.
PDRobstacle()
Construct zero-initialized.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
label size() const noexcept
Definition UList.H:706
const Cmpt & z() const noexcept
Access to the vector z component.
Definition Vector.H:145
const Cmpt & y() const noexcept
Access to the vector y component.
Definition Vector.H:140
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
Base class for surface writers.
@ BEGIN_LIST
Begin list [isseparator].
Definition token.H:174
@ END_LIST
End list [isseparator].
Definition token.H:175
@ SPACE
Space [isspace].
Definition token.H:144
An enumeration wrapper for classification of a location as being inside/outside of a volume.
Definition volumeType.H:56
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
A class for handling words, derived from Foam::string.
Definition word.H:66
static const word null
An empty word.
Definition word.H:84
OBJstream os(runTime.globalPath()/outputName)
Macros to ease declaration of member function selection tables.
#define declareMemberFunctionSelectionTable(returnType, baseType, funcName, argNames, argList, parListUnused)
Declare a run-time member-function selection (variables and adder classes).
rAUs append(new volScalarField(IOobject::groupName("rAU", phase1.name()), 1.0/(U1Eqn.A()+byDt(max(phase1.residualAlpha() - alpha1, scalar(0)) *rho1))))
Legacy and transitional routines.
labelPair readObstacleFiles(const fileName &obsFileDir, const wordList &obsFileNames, Map< obstacleGrouping > &groups)
Read obstacle files, do counting only.
Namespace for handling VTK output. Contains classes and functions for writing VTK file content.
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
Pair< label > labelPair
A pair of labels.
Definition Pair.H:54
List< word > wordList
List of word.
Definition fileName.H:60
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
MeshedSurface< face > meshedSurface
Istream & operator>>(Istream &, directionInfo &)
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
Vector< scalar > vector
Definition vector.H:57
dictionary dict
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
volScalarField & b