blockMesh.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-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::blockMesh
29 
30 Description
31  A multi-block mesh generator
32 
33  Dictionary controls
34  \table
35  Property | Description | Required | Default
36  prescale | Point scaling before transform | no | 1.0
37  scale | Point scaling after transform | no | 1.0
38  transform | Point transform (origin, rotation) | no |
39  vertices | | yes |
40  blocks | | yes |
41  edges | | no |
42  faces | | no |
43  boundary | Boundary definition | no |
44  patches | Alternate version for "boundary" | no |
45  namedBlocks | | no |
46  namedVertices | | no |
47  mergeType | Merging "points" or "topology" | no | topology
48  checkFaceCorrespondence | | no | true
49  verbose | | no | true
50  \endtable
51 
52 Note
53  The \c prescale and \c scale can be a single scalar or a vector of
54  values.
55 
56  The vertices, cells and patches for filling the blocks are demand-driven.
57 
58 SourceFiles
59  blockMesh.C
60  blockMeshCheck.C
61  blockMeshCreate.C
62  blockMeshMerge.C
63  blockMeshTopology.C
64 
65 \*---------------------------------------------------------------------------*/
66 
67 #ifndef blockMesh_H
68 #define blockMesh_H
69 
70 #include "Enum.H"
71 #include "Switch.H"
72 #include "block.H"
73 #include "PtrList.H"
74 #include "cartesianCS.H"
75 #include "searchableSurfaces.H"
76 #include "polyMesh.H"
77 #include "IOdictionary.H"
78 #include "blockVertexList.H"
79 #include "blockEdgeList.H"
80 #include "blockFaceList.H"
81 
82 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
83 
84 namespace Foam
85 {
86 
87 /*---------------------------------------------------------------------------*\
88  Class blockMesh Declaration
89 \*---------------------------------------------------------------------------*/
90 
91 class blockMesh
92 :
93  public PtrList<block>
94 {
95 public:
96 
97  // Typedefs
98 
99  //- The list of blocks is stored as a PtrList
100  typedef PtrList<block> blockList;
101 
102 
103  // Data Types
104 
105  //- The block merging strategy
106  enum mergeStrategy
107  {
108  DEFAULT_MERGE,
110  MERGE_POINTS
111  };
112 
113 
114 private:
115 
116  // Static Data
117 
118  //- Names corresponding to the mergeStrategy
119  static const Enum<mergeStrategy> strategyNames_;
120 
121 
122  // Data Types
123 
124  //- Point transformations. Internal usage
125  enum transformTypes : unsigned
126  {
127  NO_TRANSFORM = 0,
128  ROTATION = 0x1,
129  TRANSLATION = 0x2,
130  PRESCALING = 0x4,
131  PRESCALING3 = 0x8,
132  SCALING = 0x10,
133  SCALING3 = 0x20
134  };
135 
136 
137  // Private Data
138 
139  //- Reference to mesh dictionary
140  const IOdictionary& meshDict_;
141 
142  //- Output verbosity level
143  int verbose_;
144 
145  //- Check face consistency (defaults to true)
146  bool checkFaceCorrespondence_;
147 
148  //- The mesh merging strategy
149  mergeStrategy mergeStrategy_;
150 
151  //- Types of point transformations requested
152  unsigned transformType_;
153 
154  //- Optional searchable geometry to project face-points to
155  searchableSurfaces geometry_;
156 
157  //- The list of block vertices
158  blockVertexList blockVertices_;
159 
160  //- The list of block vertex positions
161  pointField vertices_;
162 
163  //- The list of curved edges
164  blockEdgeList edges_;
165 
166  //- The list of curved faces
167  blockFaceList faces_;
168 
169  //- The scaling factor before rotate/translate
170  vector prescaling_;
171 
172  //- The scaling factor after rotate/translate
173  vector scaling_;
174 
175  //- Local coordinate system transformation
176  coordSystem::cartesian transform_;
177 
178  //- The blocks themselves (the topology) as a polyMesh
179  autoPtr<polyMesh> topologyPtr_;
181  //- The sum of all cells in each block
182  label nPoints_;
184  //- The sum of all cells in each block
185  label nCells_;
186 
187  //- The merge points information
188  labelList mergeList_;
189 
190  //- The point offset added to each block. Offsets into mergeList_
191  labelList blockOffsets_;
192 
193  mutable pointField points_;
194 
195  mutable cellShapeList cells_;
196 
197  mutable faceListList patches_;
198 
199 
200  // Private Member Functions
201 
202  //- Get scaling and/or coordinate transforms
203  // \return True if scaling and/or transformations are needed
204  bool readPointTransforms(const dictionary& dict);
205 
206  void readPatches
207  (
208  const dictionary& meshDescription,
209  faceListList& tmpBlocksPatches,
212  wordList& nbrPatchNames
213  );
214 
215  void readBoundary
216  (
217  const dictionary& meshDescription,
219  faceListList& tmpBlocksPatches,
221  );
222 
223  //- Topology blocks as cell shapes
224  cellShapeList getBlockShapes() const;
225 
226  autoPtr<polyMesh> createTopology
227  (
228  const IOdictionary& dict,
229  const word& regionName
230  );
231 
232 
233  //- Simple checks for collapsed hex cells
234  bool checkDegenerate() const;
235 
236  void check(const polyMesh& bm, const dictionary& dict) const;
237 
238  //- Determine merge info and final number of cells/points
239  //- based on point distances
240  void calcGeometricalMerge();
241 
242  //- Determine merge info and final number of cells/points
243  //- based on block topology
244  void calcTopologicalMerge();
245 
246  faceList createPatchFaces(const polyPatch& patchTopologyFaces) const;
247 
248  void createPoints() const;
249  void createCells() const;
250  void createPatches() const;
251 
252 
253  //- No copy construct
254  blockMesh(const blockMesh&) = delete;
255 
256  //- No copy assignment
257  void operator=(const blockMesh&) = delete;
258 
259 
260 public:
261 
262  // Static Data
263 
264  //- The default verbosity (true)
265  static bool verboseOutput;
266 
267 
268  //- Runtime type information
269  ClassName("blockMesh");
270 
271 
272  // Constructors
273 
274  //- Construct from IOdictionary for given region
275  // Default is topological merging.
276  explicit blockMesh
277  (
278  const IOdictionary& dict,
280  mergeStrategy strategy = mergeStrategy::DEFAULT_MERGE,
281  int verbosity = 0 // 0: use static or dictionary value
282  );
283 
284 
285  //- Destructor
286  ~blockMesh() = default;
287 
288 
289  // Member Functions
290 
291  // General Access, Description
292 
293  //- Access to input dictionary
294  const dictionary& meshDict() const noexcept
295  {
296  return meshDict_;
297  }
298 
299  //- Optional searchable geometry to project face-points to
300  const searchableSurfaces& geometry() const noexcept
301  {
302  return geometry_;
303  }
304 
305  //- The curved edges
306  const blockEdgeList& edges() const noexcept
307  {
308  return edges_;
309  }
310 
311  //- The curved faces
312  const blockFaceList& faces() const noexcept
313  {
314  return faces_;
315  }
316 
317  //- True if the blockMesh topology exists
318  bool good() const noexcept;
319 
320  //- Same as good()
321  bool valid() const noexcept { return good(); }
322 
323  //- Return the patch names
324  wordList patchNames() const;
325 
326  //- Number of blocks with specified zones
327  label numZonedBlocks() const;
328 
329 
330  // Point Transformations
331 
332  //- True if scaling and/or transformations are needed
333  bool hasPointTransforms() const noexcept;
334 
335  //- Apply coordinate transforms and scaling
337 
338  //- Apply coordinate transforms and scaling
339  tmp<pointField> globalPosition(const pointField& localPoints) const;
340 
341 
342  // Block Topology
343 
344  //- Reference to point field defining the blocks,
345  //- these points are \b unscaled and \b non-transformed
346  const pointField& vertices() const noexcept;
347 
348  //- Point field defining the blocks,
349  //- optionally transformed and scaled
350  tmp<pointField> vertices(bool applyTransform) const;
351 
352  //- The blockMesh topology as a polyMesh
353  //- \b unscaled and \b non-transformed
354  const polyMesh& topology() const;
355 
356  //- The blockMesh topology as a polyMesh
357  //- optionally transformed and scaled
358  refPtr<polyMesh> topology(bool applyTransform) const;
359 
360 
361  // Detailed Mesh
362 
363  //- The points for the entire mesh.
364  //- These points \b are scaled and transformed
365  const pointField& points() const;
366 
367  //- Return cell shapes list
368  const cellShapeList& cells() const;
369 
370  //- Return the patch face lists
371  const faceListList& patches() const;
372 
373  //- Patch information from the topology mesh
374  PtrList<dictionary> patchDicts() const;
375 
376 
377  // Verbosity
378 
379  //- Output verbosity level
380  int verbose() const noexcept;
381 
382  //- Change the output verbosity level.
383  // \return old level
384  int verbose(const int level) noexcept;
385 
386 
387  // Mesh Generation
388 
389  //- Create polyMesh, with cell zones
390  autoPtr<polyMesh> mesh(const IOobject& io) const;
391 
392 
393  // Housekeeping
394 
395  //- Old (v2106 and earlier) uniform scaling factor
396  // \deprecated use inplacePointTransforms or globalPosition instead
397  scalar scaleFactor() const
398  {
399  return scaling_.x();
400  }
401 };
402 
403 
404 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
405 
406 } // End namespace Foam
407 
408 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
409 
410 #endif
411 
412 // ************************************************************************* //
wordList patchNames() const
Return the patch names.
Definition: blockMesh.C:392
dictionary dict
List< faceList > faceListList
List of faceList.
Definition: faceListFwd.H:44
~blockMesh()=default
Destructor.
const pointField & vertices() const noexcept
Reference to point field defining the blocks, these points are unscaled and non-transformed.
Definition: blockMesh.C:320
static bool verboseOutput
The default verbosity (true)
Definition: blockMesh.H:396
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const blockEdgeList & edges() const noexcept
The curved edges.
Definition: blockMesh.H:450
PtrList< blockFace > blockFaceList
A PtrList of blockFaces.
Definition: blockFaceList.H:41
wordList patchTypes(nPatches)
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:381
List< cellShape > cellShapeList
List of cellShape.
Definition: cellShapeList.H:32
scalar scaleFactor() const
Old (v2106 and earlier) uniform scaling factor.
Definition: blockMesh.H:583
bool good() const noexcept
True if the blockMesh topology exists.
Definition: blockMesh.C:300
tmp< pointField > globalPosition(const pointField &localPoints) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:507
"topology" merge by block topology (default)
Definition: blockMesh.H:183
const searchableSurfaces & geometry() const noexcept
Optional searchable geometry to project face-points to.
Definition: blockMesh.H:442
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:50
A Cartesian coordinate system.
Definition: cartesianCS.H:65
PtrList< block > blockList
The list of blocks is stored as a PtrList.
Definition: blockMesh.H:172
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
A multi-block mesh generator.
Definition: blockMesh.H:161
mergeStrategy
The block merging strategy.
Definition: blockMesh.H:180
A class for handling words, derived from Foam::string.
Definition: word.H:63
const pointField & points() const
The points for the entire mesh. These points are scaled and transformed.
Definition: blockMesh.C:359
constexpr PtrList() noexcept
Default construct.
Definition: PtrListI.H:29
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
Default (TOPOLOGY), not selectable.
Definition: blockMesh.H:182
Container for searchableSurfaces. The collection is specified as a dictionary. For example...
const polyMesh & topology() const
The blockMesh topology as a polyMesh unscaled and non-transformed.
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
"points" merge by point geometry
Definition: blockMesh.H:184
PtrList< blockVertex > blockVertexList
A PtrList of blockVertex.
const dictionary & meshDict() const noexcept
Access to input dictionary.
Definition: blockMesh.H:434
List< word > wordList
List of word.
Definition: fileName.H:59
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:370
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
PtrList< blockEdge > blockEdgeList
A PtrList of blockEdges.
Definition: blockEdgeList.H:41
bool valid() const noexcept
Same as good()
Definition: blockMesh.H:471
bool hasPointTransforms() const noexcept
True if scaling and/or transformations are needed.
Definition: blockMesh.C:428
autoPtr< polyMesh > mesh(const IOobject &io) const
Create polyMesh, with cell zones.
bool inplacePointTransforms(pointField &pts) const
Apply coordinate transforms and scaling.
Definition: blockMesh.C:434
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
PtrList< dictionary > patchDicts() const
Patch information from the topology mesh.
Definition: blockMesh.C:342
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
ClassName("blockMesh")
Runtime type information.
int verbose() const noexcept
Output verbosity level.
Definition: blockMesh.C:306
Namespace for OpenFOAM.
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:410
const pointField & pts
const blockFaceList & faces() const noexcept
The curved faces.
Definition: blockMesh.H:458