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 -------------------------------------------------------------------------------
11 License
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 
27 Class
28  Foam::PDRobstacle
29 
30 Description
31  Obstacle definitions for PDR
32 
33 SourceFiles
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 
55 namespace Foam
56 {
57 
58 // Forward Declarations
59 class boundBox;
60 class PDRobstacle;
61 
62 Istream& operator>>(Istream& is, PDRobstacle& obs);
63 Ostream& operator<<(Ostream& os, const InfoProxy<PDRobstacle>& info);
64 
65 namespace vtk
66 {
67  class surfaceWriter;
68 }
69 
70 
71 /*---------------------------------------------------------------------------*\
72  Class PDRobstacle Declaration
73 \*---------------------------------------------------------------------------*/
74 
75 class PDRobstacle
76 {
77 public:
78 
79  //- Obstacle types (legacy numbering)
81  {
82  NONE = 0,
83  CUBOID_1 = 1,
84  CYLINDER = 2,
87  CUBOID = 6,
88  WALL_BEAM = 7,
89  GRATING = 8,
90  OLD_INLET = 9,
91  OLD_BLOWOFF = 10,
92  CIRC_PATCH = 12,
93  RECT_PATCH = 16,
94  DIAG_BEAM = 22,
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
140  {
141  scalar wa;
142  scalar slat_width;
143  scalar blowoff_press;
144  };
145  union
146  {
147  scalar wb;
148  scalar blowoff_time;
149  };
150  scalar vbkge;
151  scalar xbkge;
152  scalar ybkge;
153  scalar zbkge;
155  union
156  {
157  int blowoff_type;
158  int inlet_dirn;
159  };
161  string identifier;
162 
163 public:
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,
179  PDRobstacle,
180  read,
181  dictionary,
182  (
183  PDRobstacle& obs,
184  const dictionary& dict
185  ),
186  (obs, dict)
187  );
188 
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
200  DynamicList<PDRobstacle>& blocks,
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
212  DynamicList<PDRobstacle>& blocks,
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)
259  (
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
276  (
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 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 // Global Operators
307 
308 //- Compare according to x0 position
309 bool operator<(const PDRobstacle& a, const PDRobstacle& b);
310 
311 //- For list output, assert that no obstacles are identical
312 inline bool operator!=(const PDRobstacle& a, const PDRobstacle& b)
313 {
314  return true;
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 namespace PDRlegacy
321 {
322 
323 /*---------------------------------------------------------------------------*\
324  Class obstacleGrouping Declaration
325 \*---------------------------------------------------------------------------*/
326 
327 //- Locations for each instance of an obstacle group.
328 class 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_;
337 
338 
339 public:
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()
372  {
373  ++nCylinder_;
374  }
375 
376  //- The number of obstacles
377  label nObstacle() const
378  {
379  return nObstacle_;
380  }
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
390  {
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,
425  Map<obstacleGrouping>& groups
426 );
427 
428 
429 //- Read obstacle files and add to the lists
430 // \return the total volume
431 scalar readObstacleFiles
432 (
433  const fileName& obsFileDir,
434  const wordList& obsFileNames,
435  const Map<obstacleGrouping>& groups,
436  const boundBox& meshBb,
437  // output
438  DynamicList<PDRobstacle>& blocks,
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.
453 inline Ostream& operator<<
454 (
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 // ************************************************************************* //
dictionary dict
scalar volume() const
Volume of the obstacle.
void addObstacle()
Increment the number of obstacles.
Definition: PDRobstacle.H:439
uint8_t direction
Definition: direction.H:46
A class for handling file names.
Definition: fileName.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label nCylinder() const
The number of cylinder-like obstacles.
Definition: PDRobstacle.H:463
labelPair readObstacleFiles(const fileName &obsFileDir, const wordList &obsFileNames, Map< obstacleGrouping > &groups)
Read obstacle files, do counting only.
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.
An Istream is an abstract base class for all input systems (streams, files, token lists etc)...
Definition: Istream.H:57
Write faces/points (optionally with fields) as a vtp file or a legacy vtk file.
scalar sortBias
Bias for position sorting.
Definition: PDRobstacle.H:128
static label addPieces(vtk::surfaceWriter &surfWriter, const UList< PDRobstacle > &list, label pieceId=0)
Add pieces to vtp output.
Begin list [isseparator].
Definition: token.H:161
void clear()
Reset to a zero obstacle.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
static constexpr int maxBlowoffPressure
The max blowoff pressure [bar].
Definition: PDRobstacle.H:105
An enumeration wrapper for classification of a location as being inside/outside of a volume...
Definition: volumeType.H:55
bool operator<(const IOstreamOption::versionNumber &a, const IOstreamOption::versionNumber &b) noexcept
Version A older than B.
label nObstacle() const
The number of obstacles.
Definition: PDRobstacle.H:455
declareMemberFunctionSelectionTable(void, PDRobstacle, read, dictionary,(PDRobstacle &obs, const dictionary &dict),(obs, dict))
constexpr const char *const group
Group name for atomic constants.
void append(const scalar x, const scalar y, const scalar z)
Add a location.
Definition: PDRobstacle.H:500
scalar y
InfoProxy< PDRobstacle > info() const noexcept
Return info proxy, used to print information to a stream.
Definition: PDRobstacle.H:351
scalar z() const
Definition: PDRobstacle.H:257
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:51
scalar y() const
Definition: PDRobstacle.H:256
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
scalar dia() const
Definition: PDRobstacle.H:144
Istream & operator>>(Istream &, directionInfo &)
scalar len() const
Definition: PDRobstacle.H:146
Locations for each instance of an obstacle group.
Definition: PDRobstacle.H:390
Space [isspace].
Definition: token.H:131
static const word null
An empty word.
Definition: word.H:84
label groupId
The group-id.
Definition: PDRobstacle.H:113
PDRobstacle()
Construct zero-initialized.
End list [isseparator].
Definition: token.H:162
obstacleGrouping()
Construct null.
Definition: PDRobstacle.H:410
Vector< scalar > vector
Definition: vector.H:57
void append(const point &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
label nTotalCylinder() const
The number of locations x number of cylinder-like obstacles.
Definition: PDRobstacle.H:479
bool isCylinder() const
Is obstacle cylinder-like?
Definition: PDRobstacleI.H:33
volumeType trim(const boundBox &bb)
Trim obstacle to ensure it is within the specified bounding box and return the intersection type...
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:105
void addCylinder()
Increment the number of cylinder-like obstacles.
Definition: PDRobstacle.H:447
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
const direction noexcept
Definition: Scalar.H:258
static void generateVtk(const fileName &outputDir, const UList< PDRobstacle > &obslist, const UList< PDRobstacle > &cyllist)
Generate multi-piece VTK (vtp) file of obstacles.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
label nTotalObstacle() const
The number of locations x number of obstacles.
Definition: PDRobstacle.H:471
OBJstream os(runTime.globalPath()/outputName)
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
Definition: DynamicListI.H:405
label nTotal() const
The number of locations x number of obstacles.
Definition: PDRobstacle.H:487
void readProperties(const dictionary &dict)
Read the &#39;name&#39; identifier if present.
Macros to ease declaration of member function selection tables.
A helper class for outputting values to Ostream.
Definition: ensightCells.H:43
vector point
Point is a vector.
Definition: point.H:37
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.
legacyTypes
Obstacle types (legacy numbering)
Definition: PDRobstacle.H:77
PtrList< volScalarField > & Y
vector span
The obstacle dimensions (for boxes)
Definition: PDRobstacle.H:140
int typeId
The obstacle type-id.
Definition: PDRobstacle.H:118
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 operator!=(const eddy &a, const eddy &b)
Definition: eddy.H:297
bool tooSmall(const scalar minWidth) const
True if the obstacle is considered to be too small.
direction orient
The x/y/z orientation (0,1,2)
Definition: PDRobstacle.H:123
void scale(const scalar factor)
Scale obstacle dimensions by specified scaling factor.
friend Istream & operator>>(Istream &is, PDRobstacle &obs)
point pt
The obstacle location.
Definition: PDRobstacle.H:135
scalar x() const
Obstacle position accessors.
Definition: PDRobstacle.H:255
Obstacle definitions for PDR.
Definition: PDRobstacle.H:70
meshedSurface surface() const
Surface (points, faces) representation.
List< treeBoundBox > meshBb(1, treeBoundBox(coarseMesh.points()).extend(rndGen, 1e-3))
Namespace for OpenFOAM.
void clear()
Clear obstacle count and locations.
Definition: PDRobstacle.H:429
scalar theta() const
Definition: PDRobstacle.H:145
A HashTable to objects of type <T> with a label key.
bool read(Istream &is)
Read name / dictionary.