blockMesh.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::blockMesh
26 
27 Description
28  A multi-block mesh generator
29 
30 Note
31  The vertices, cells and patches for filling the blocks are demand-driven.
32 
33 SourceFiles
34  blockMesh.C
35  blockMeshCheck.C
36  blockMeshCreate.C
37  blockMeshMerge.C
38  blockMeshTopology.C
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #ifndef blockMesh_H
43 #define blockMesh_H
44 
45 #include "blockList.H"
46 #include "searchableSurfaces.H"
47 #include "polyMesh.H"
48 #include "IOdictionary.H"
49 #include "blockVertexList.H"
50 #include "blockEdgeList.H"
51 #include "blockFaceList.H"
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 namespace Foam
56 {
57 
58 /*---------------------------------------------------------------------------*\
59  Class blockMesh Declaration
60 \*---------------------------------------------------------------------------*/
61 
62 class blockMesh
63 :
64  public blockList
65 {
66  // Private data
67 
68  //- Reference to mesh dictionary
69  const IOdictionary& meshDict_;
70 
71  //- Switch for verbose output
72  bool verboseOutput;
73 
74  //- Switch checking face consistency (defaults to true)
75  bool checkFaceCorrespondence_;
76 
77  //- Optional searchable geometry to project face-points to
78  searchableSurfaces geometry_;
79 
80  //- The scaling factor to convert to meters
81  scalar scaleFactor_;
82 
83  //- The list of block vertices
84  blockVertexList blockVertices_;
85 
86  //- The list of block vertex positions
87  pointField vertices_;
88 
89  //- The list of curved edges
90  blockEdgeList edges_;
91 
92  //- The list of curved faces
93  blockFaceList faces_;
94 
95  //- The blocks themselves (the topology) as a polyMesh
96  polyMesh* topologyPtr_;
97 
98  label nPoints_;
99 
100  //- The sum of all cells in each block
101  label nCells_;
102 
103  //- The point offset added to each block
104  labelList blockOffsets_;
105 
106  //- The merge points information
107  labelList mergeList_;
108 
109  mutable pointField points_;
110 
111  mutable cellShapeList cells_;
112 
113  mutable faceListList patches_;
114 
115 
116  // Private Member Functions
117 
118  template<class Source>
119  void checkPatchLabels
120  (
121  const Source& source,
122  const word& patchName,
123  const pointField& points,
124  faceList& patchShapes
125  ) const;
126 
127  void readPatches
128  (
129  const dictionary& meshDescription,
130  faceListList& tmpBlocksPatches,
133  wordList& nbrPatchNames
134  );
135 
136  void readBoundary
137  (
138  const dictionary& meshDescription,
139  wordList& patchNames,
140  faceListList& tmpBlocksPatches,
142  );
143 
144  void createCellShapes(cellShapeList& tmpBlockCells);
145 
146  polyMesh* createTopology(const IOdictionary&, const word& regionName);
147 
148  void check(const polyMesh&, const dictionary&) const;
149 
150  //- Determine the merge info and the final number of cells/points
151  void calcMergeInfo();
152 
153  //- Determine the merge info and the final number of cells/points
154  void calcMergeInfoFast();
155 
156  faceList createPatchFaces(const polyPatch& patchTopologyFaces) const;
157 
158  void createPoints() const;
159  void createCells() const;
160  void createPatches() const;
161 
162  //- As copy (not implemented)
163  blockMesh(const blockMesh&);
164 
165 
166 public:
167 
168  // Static data members
169 
170  ClassName("blockMesh");
171 
172 
173  // Constructors
174 
175  //- Construct from IOdictionary
176  blockMesh(const IOdictionary&, const word& regionName);
177 
178 
179  //- Destructor
180  ~blockMesh();
181 
182 
183  // Member Functions
184 
185  // Access
186 
187  //- Access to input dictionary
188  const dictionary& meshDict() const
189  {
190  return meshDict_;
191  }
192 
193  //- Optional searchable geometry to project face-points to
194  const searchableSurfaces& geometry() const
195  {
196  return geometry_;
197  }
198 
199  //- Reference to point field defining the blockMesh
200  // these points have not been scaled by scaleFactor
201  const pointField& vertices() const;
202 
203  //- Return the blockMesh topology as a polyMesh
204  const polyMesh& topology() const;
205 
206  //- Return the curved edges
207  const blockEdgeList& edges() const
208  {
209  return edges_;
210  }
211 
212  //- Return the curved faces
213  const blockFaceList& faces() const
214  {
215  return faces_;
216  }
217 
218  //- The scaling factor used to convert to meters
219  scalar scaleFactor() const;
220 
221  //- The points for the entire mesh
222  // these points have been scaled by scaleFactor
223  const pointField& points() const;
224 
225  //- Return cell shapes list
226  const cellShapeList& cells() const;
227 
228  //- Return the patch face lists
229  const faceListList& patches() const;
230 
231  //- Get patch information from the topology mesh
233 
234  //- Return patch names
235  wordList patchNames() const;
236 
237  //- Number of blocks with specified zones
238  label numZonedBlocks() const;
239 
240 
241  // Edit
242 
243  //- Enable/disable verbose information about the progress
244  void verbose(const bool on=true);
245 
246 
247  // Write
248 
249  //- Writes edges of blockMesh in OBJ format.
250  void writeTopology(Ostream&) const;
251 };
252 
253 
254 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
255 
256 } // End namespace Foam
257 
258 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
259 
260 #endif
261 
262 // ************************************************************************* //
const blockFaceList & faces() const
Return the curved faces.
Definition: blockMesh.H:212
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Foam::word regionName
const faceListList & patches() const
Return the patch face lists.
wordList patchTypes(nPatches)
void writeTopology(Ostream &) const
Writes edges of blockMesh in OBJ format.
scalar scaleFactor() const
The scaling factor used to convert to meters.
const pointField & vertices() const
Reference to point field defining the blockMesh.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
A multi-block mesh generator.
Definition: blockMesh.H:61
A class for handling words, derived from string.
Definition: word.H:59
const searchableSurfaces & geometry() const
Optional searchable geometry to project face-points to.
Definition: blockMesh.H:193
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Container for searchableSurfaces.
const blockEdgeList & edges() const
Return the curved edges.
Definition: blockMesh.H:206
void verbose(const bool on=true)
Enable/disable verbose information about the progress.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
PtrList< dictionary > patchDicts() const
Get patch information from the topology mesh.
~blockMesh()
Destructor.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
const cellShapeList & cells() const
Return cell shapes list.
const pointField & points() const
The points for the entire mesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label numZonedBlocks() const
Number of blocks with specified zones.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const dictionary & meshDict() const
Access to input dictionary.
Definition: blockMesh.H:187
wordList patchNames() const
Return patch names.
ClassName("blockMesh")
Namespace for OpenFOAM.