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 | Copyright (C) 2011-2016 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  //- Optional searchable geometry to project face-points to
75  searchableSurfaces geometry_;
76 
77  //- The scaling factor to convert to metres
78  scalar scaleFactor_;
79 
80  //- The list of block vertices
81  blockVertexList blockVertices_;
82 
83  //- The list of block vertex positions
84  pointField vertices_;
85 
86  //- The list of curved edges
87  blockEdgeList edges_;
88 
89  //- The list of curved faces
90  blockFaceList faces_;
91 
92  //- The blocks themselves (the topology) as a polyMesh
93  polyMesh* topologyPtr_;
94 
95  label nPoints_;
96 
97  //- The sum of all cells in each block
98  label nCells_;
99 
100  //- The point offset added to each block
101  labelList blockOffsets_;
102 
103  //- The merge points information
104  labelList mergeList_;
105 
106  mutable pointField points_;
107 
108  mutable cellShapeList cells_;
109 
110  mutable faceListList patches_;
111 
112 
113  // Private Member Functions
114 
115  template<class Source>
116  void checkPatchLabels
117  (
118  const Source& source,
119  const word& patchName,
120  const pointField& points,
121  faceList& patchShapes
122  ) const;
123 
124  void readPatches
125  (
126  const dictionary& meshDescription,
127  faceListList& tmpBlocksPatches,
130  wordList& nbrPatchNames
131  );
132 
133  void readBoundary
134  (
135  const dictionary& meshDescription,
136  wordList& patchNames,
137  faceListList& tmpBlocksPatches,
139  );
140 
141  void createCellShapes(cellShapeList& tmpBlockCells);
142 
143  polyMesh* createTopology(const IOdictionary&, const word& regionName);
144 
145  void check(const polyMesh&, const dictionary&) const;
146 
147  //- Determine the merge info and the final number of cells/points
148  void calcMergeInfo();
149 
150  //- Determine the merge info and the final number of cells/points
151  void calcMergeInfoFast();
152 
153  faceList createPatchFaces(const polyPatch& patchTopologyFaces) const;
154 
155  void createPoints() const;
156  void createCells() const;
157  void createPatches() const;
158 
159  //- As copy (not implemented)
160  blockMesh(const blockMesh&);
161 
162 
163 public:
164 
165  // Static data members
166 
167  ClassName("blockMesh");
168 
169 
170  // Constructors
171 
172  //- Construct from IOdictionary
173  blockMesh(const IOdictionary&, const word& regionName);
174 
175 
176  //- Destructor
177  ~blockMesh();
178 
179 
180  // Member Functions
181 
182  // Access
183 
184  //- Access to input dictionary
185  const dictionary& meshDict() const
186  {
187  return meshDict_;
188  }
189 
190  //- Optional searchable geometry to project face-points to
191  const searchableSurfaces& geometry() const
192  {
193  return geometry_;
194  }
195 
196  //- Reference to point field defining the blockMesh
197  // these points have not been scaled by scaleFactor
198  const pointField& vertices() const;
199 
200  //- Return the blockMesh topology as a polyMesh
201  const polyMesh& topology() const;
202 
203  //- Return the curved edges
204  const blockEdgeList& edges() const
205  {
206  return edges_;
207  }
208 
209  //- Return the curved faces
210  const blockFaceList& faces() const
211  {
212  return faces_;
213  }
214 
215  //- The scaling factor used to convert to metres
216  scalar scaleFactor() const;
217 
218  //- The points for the entire mesh
219  // these points have been scaled by scaleFactor
220  const pointField& points() const;
221 
222  //- Return cell shapes list
223  const cellShapeList& cells() const;
224 
225  //- Return the patch face lists
226  const faceListList& patches() const;
227 
228  //- Get patch information from the topology mesh
230 
231  //- Return patch names
232  wordList patchNames() const;
233 
234  //- Number of blocks with specified zones
235  label numZonedBlocks() const;
236 
237 
238  // Edit
239 
240  //- Enable/disable verbose information about the progress
241  void verbose(const bool on=true);
242 
243 
244  // Write
245 
246  //- Writes edges of blockMesh in OBJ format.
247  void writeTopology(Ostream&) const;
248 };
249 
250 
251 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
252 
253 } // End namespace Foam
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 #endif
258 
259 // ************************************************************************* //
wordList patchNames() const
Return patch names.
Definition: blockMesh.C:172
const blockFaceList & faces() const
Return the curved faces.
Definition: blockMesh.H:209
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
wordList patchTypes(nPatches)
void writeTopology(Ostream &) const
Writes edges of blockMesh in OBJ format.
Definition: blockMesh.C:206
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:161
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
scalar scaleFactor() const
The scaling factor used to convert to metres.
Definition: blockMesh.C:133
A multi-block mesh generator.
Definition: blockMesh.H:61
A class for handling words, derived from string.
Definition: word.H:59
const pointField & points() const
The points for the entire mesh.
Definition: blockMesh.C:139
const searchableSurfaces & geometry() const
Optional searchable geometry to project face-points to.
Definition: blockMesh.H:190
Container for searchableSurfaces.
const blockEdgeList & edges() const
Return the curved edges.
Definition: blockMesh.H:203
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Definition: blockMesh.C:103
void verbose(const bool on=true)
Enable/disable verbose information about the progress.
Definition: blockMesh.C:91
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
~blockMesh()
Destructor.
Definition: blockMesh.C:83
const pointField & vertices() const
Reference to point field defining the blockMesh.
Definition: blockMesh.C:97
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:150
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
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
PtrList< dictionary > patchDicts() const
Get patch information from the topology mesh.
Definition: blockMesh.C:116
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:184
ClassName("blockMesh")
Namespace for OpenFOAM.
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:190