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-2015 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 "polyMesh.H"
47 #include "IOdictionary.H"
48 #include "curvedEdgeList.H"
49 
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
51 
52 namespace Foam
53 {
54 
55 /*---------------------------------------------------------------------------*\
56  Class blockMesh Declaration
57 \*---------------------------------------------------------------------------*/
58 
59 class blockMesh
60 :
61  public blockList
62 {
63  // Private data
64 
65  //- Switch for verbose output
66  static bool verboseOutput;
67 
68  //- Point field defining the block mesh (corners)
70 
71  //- The list of curved edges
73 
74  //- The scaling factor to convert to metres
75  scalar scaleFactor_;
76 
77  //- The blocks themselves (the topology) as a polyMesh
79 
81 
82  //- The sum of all cells in each block
83  label nCells_;
84 
85  //- The point offset added to each block
87 
88  //- The merge points information
90 
91  mutable pointField points_;
92 
93  mutable cellShapeList cells_;
94 
95  mutable faceListList patches_;
96 
97 
98  // Private Member Functions
99 
100  bool blockLabelsOK
101  (
102  const label blockLabel,
103  const pointField& points,
104  const cellShape& blockShape
105  ) const;
106 
107  bool patchLabelsOK
108  (
109  const label patchLabel,
110  const pointField& points,
111  const faceList& patchShapes
112  ) const;
113 
114  bool readPatches
115  (
116  const dictionary& meshDescription,
117  faceListList& tmpBlocksPatches,
120  wordList& nbrPatchNames
121  );
122 
123  bool readBoundary
124  (
125  const dictionary& meshDescription,
127  faceListList& tmpBlocksPatches,
129  );
130 
131  void createCellShapes(cellShapeList& tmpBlockCells);
132 
134  void checkBlockMesh(const polyMesh&) const;
135 
136  //- Determine the merge info and the final number of cells/points
137  void calcMergeInfo();
138 
139  //- Determine the merge info and the final number of cells/points
140  void calcMergeInfoFast();
141 
142  faceList createPatchFaces(const polyPatch& patchTopologyFaces) const;
143 
144  void createPoints() const;
145  void createCells() const;
146  void createPatches() const;
147 
148  //- As copy (not implemented)
149  blockMesh(const blockMesh&);
150 
151 
152 public:
153 
154  // Static data members
155 
156  ClassName("blockMesh");
157 
158 
159  // Constructors
160 
161  //- Construct from IOdictionary
162  blockMesh(const IOdictionary&, const word& regionName);
163 
164 
165  //- Destructor
166  ~blockMesh();
167 
168 
169  // Member Functions
170 
171  // Access
172 
173  //- Reference to point field defining the block mesh
174  // these points have not been scaled by scaleFactor
175  const pointField& blockPointField() const;
176 
177  //- Return the blockMesh topology as a polyMesh
178  const polyMesh& topology() const;
179 
180  //- Return the curved edges
181  const curvedEdgeList& edges() const
182  {
183  return edges_;
184  }
185 
186  //- The scaling factor used to convert to metres
187  scalar scaleFactor() const;
188 
189  //- The points for the entire mesh
190  // these points have been scaled by scaleFactor
191  const pointField& points() const;
192 
193  //- Return cell shapes list
194  const cellShapeList& cells() const;
195 
196  //- Return the patch face lists
197  const faceListList& patches() const;
198 
199  //- Get patch information from the topology mesh
201 
202  //- Return patch names
203  wordList patchNames() const;
204 
205  //- Number of blocks with specified zones
206  label numZonedBlocks() const;
207 
208 
209  // Edit
210 
211  //- Clear geometry (internal points, cells, boundaryPatches)
212  void clearGeom();
213 
214  //- Enable/disable verbose information about the progress
215  static void verbose(const bool on=true);
216 
217 
218  // Write
219 
220  //- Writes edges of blockMesh in OBJ format.
221  void writeTopology(Ostream&) const;
222 };
223 
224 
225 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
226 
227 } // End namespace Foam
228 
229 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
230 
231 #endif
232 
233 // ************************************************************************* //
curvedEdgeList.H
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::blockMesh::scaleFactor_
scalar scaleFactor_
The scaling factor to convert to metres.
Definition: blockMesh.H:74
Foam::blockMesh::patchDicts
PtrList< dictionary > patchDicts() const
Get patch information from the topology mesh.
Definition: blockMesh.C:95
Foam::blockMesh::readPatches
bool readPatches(const dictionary &meshDescription, faceListList &tmpBlocksPatches, wordList &patchNames, wordList &patchTypes, wordList &nbrPatchNames)
Definition: blockMeshTopology.C:34
Foam::blockMesh::patches
const faceListList & patches() const
Return the patch face lists.
Definition: blockMesh.C:140
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::blockMesh::ClassName
ClassName("blockMesh")
Foam::blockMesh::readBoundary
bool readBoundary(const dictionary &meshDescription, wordList &patchNames, faceListList &tmpBlocksPatches, PtrList< dictionary > &patchDicts)
Definition: blockMeshTopology.C:186
Foam::blockMesh::checkBlockMesh
void checkBlockMesh(const polyMesh &) const
Definition: blockMeshCheck.C:30
Foam::blockMesh::blockPointField_
pointField blockPointField_
Point field defining the block mesh (corners)
Definition: blockMesh.H:68
Foam::blockMesh::nPoints_
label nPoints_
Definition: blockMesh.H:79
polyMesh.H
Foam::blockMesh::topology
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Definition: blockMesh.C:82
Foam::blockMesh::patches_
faceListList patches_
Definition: blockMesh.H:94
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::blockMesh::verbose
static void verbose(const bool on=true)
Enable/disable verbose information about the progress.
Definition: blockMesh.C:70
patchTypes
wordList patchTypes(nPatches)
Foam::blockMesh::createCellShapes
void createCellShapes(cellShapeList &tmpBlockCells)
Definition: blockMeshTopology.C:235
Foam::blockMesh::writeTopology
void writeTopology(Ostream &) const
Writes edges of blockMesh in OBJ format.
Definition: blockMesh.C:185
Foam::blockMesh::blockLabelsOK
bool blockLabelsOK(const label blockLabel, const pointField &points, const cellShape &blockShape) const
Definition: blockMeshCheck.C:150
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::blockMesh::nCells_
label nCells_
The sum of all cells in each block.
Definition: blockMesh.H:82
Foam::blockMesh::calcMergeInfoFast
void calcMergeInfoFast()
Determine the merge info and the final number of cells/points.
Definition: blockMeshMergeFast.C:291
Foam::label
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
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::blockMesh::createCells
void createCells() const
Definition: blockMeshCreate.C:99
Foam::blockMesh::~blockMesh
~blockMesh()
Destructor.
Definition: blockMesh.C:62
Foam::blockMesh::edges
const curvedEdgeList & edges() const
Return the curved edges.
Definition: blockMesh.H:180
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::blockMesh::edges_
curvedEdgeList edges_
The list of curved edges.
Definition: blockMesh.H:71
Foam::blockMesh::clearGeom
void clearGeom()
Clear geometry (internal points, cells, boundaryPatches)
Definition: blockMeshCreate.C:283
Foam::blockMesh::scaleFactor
scalar scaleFactor() const
The scaling factor used to convert to metres.
Definition: blockMesh.C:112
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::blockMesh::createPatches
void createPatches() const
Definition: blockMeshCreate.C:259
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::blockMesh::points_
pointField points_
Definition: blockMesh.H:90
Foam::blockMesh::createPoints
void createPoints() const
Definition: blockMeshCreate.C:32
IOdictionary.H
Foam::blockMesh::numZonedBlocks
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:169
Foam::blockMesh::points
const pointField & points() const
The points for the entire mesh.
Definition: blockMesh.C:118
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::blockMesh::patchNames
wordList patchNames() const
Return patch names.
Definition: blockMesh.C:151
Foam::blockMesh::blockPointField
const pointField & blockPointField() const
Reference to point field defining the block mesh.
Definition: blockMesh.C:76
Foam::blockMesh::blockMesh
blockMesh(const blockMesh &)
As copy (not implemented)
Foam::blockMesh::blockOffsets_
labelList blockOffsets_
The point offset added to each block.
Definition: blockMesh.H:85
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::blockMesh::mergeList_
labelList mergeList_
The merge points information.
Definition: blockMesh.H:88
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:129
Foam::blockMesh::createPatchFaces
faceList createPatchFaces(const polyPatch &patchTopologyFaces) const
Definition: blockMeshCreate.C:145
Foam::blockMesh::createTopology
polyMesh * createTopology(const IOdictionary &, const word &regionName)
Definition: blockMeshTopology.C:259
Foam::blockMesh::verboseOutput
static bool verboseOutput
Switch for verbose output.
Definition: blockMesh.H:65
blockList.H
Foam::blockMesh::calcMergeInfo
void calcMergeInfo()
Determine the merge info and the final number of cells/points.
Definition: blockMeshMerge.C:30
Foam::blockMesh::cells_
cellShapeList cells_
Definition: blockMesh.H:92
Foam::blockMesh::topologyPtr_
polyMesh * topologyPtr_
The blocks themselves (the topology) as a polyMesh.
Definition: blockMesh.H:77
Foam::blockMesh::patchLabelsOK
bool patchLabelsOK(const label patchLabel, const pointField &points, const faceList &patchShapes) const
Definition: blockMeshCheck.C:185
Foam::blockMesh
A multi-block mesh generator.
Definition: blockMesh.H:58