blockMesh.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "blockMesh.H"
27 #include "Time.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33  defineDebugSwitch(blockMesh, 0);
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
38 
39 Foam::blockMesh::blockMesh(const IOdictionary& dict, const word& regionName)
40 :
41  meshDict_(dict),
42  verboseOutput(meshDict_.lookupOrDefault<Switch>("verbose", true)),
43  geometry_
44  (
45  IOobject
46  (
47  "geometry", // dummy name
48  meshDict_.time().constant(), // instance
49  "geometry", // local
50  meshDict_.time(), // registry
51  IOobject::MUST_READ,
52  IOobject::NO_WRITE
53  ),
54  meshDict_.found("geometry")
55  ? meshDict_.subDict("geometry")
56  : dictionary(),
57  true
58  ),
59  scaleFactor_(1.0),
60  blockVertices_
61  (
62  meshDict_.lookup("vertices"),
63  blockVertex::iNew(meshDict_, geometry_)
64  ),
65  vertices_(Foam::vertices(blockVertices_)),
66  topologyPtr_(createTopology(meshDict_, regionName))
67 {
68  Switch fastMerge(meshDict_.lookupOrDefault<Switch>("fastMerge", false));
69 
70  if (fastMerge)
71  {
72  calcMergeInfoFast();
73  }
74  else
75  {
76  calcMergeInfo();
77  }
78 }
79 
80 
81 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
82 
84 {
85  delete topologyPtr_;
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
91 void Foam::blockMesh::verbose(const bool on)
92 {
93  verboseOutput = on;
94 }
95 
96 
98 {
99  return vertices_;
100 }
101 
102 
104 {
105  if (!topologyPtr_)
106  {
108  << "topologyPtr_ not allocated"
109  << exit(FatalError);
110  }
111 
112  return *topologyPtr_;
113 }
114 
115 
117 {
118  const polyPatchList& patchTopologies = topology().boundaryMesh();
119 
120  PtrList<dictionary> patchDicts(patchTopologies.size());
121 
122  forAll(patchTopologies, patchi)
123  {
124  OStringStream os;
125  patchTopologies[patchi].write(os);
126  IStringStream is(os.str());
127  patchDicts.set(patchi, new dictionary(is));
128  }
129  return patchDicts;
130 }
131 
132 
133 Foam::scalar Foam::blockMesh::scaleFactor() const
134 {
135  return scaleFactor_;
136 }
137 
138 
140 {
141  if (points_.empty())
142  {
143  createPoints();
144  }
145 
146  return points_;
147 }
148 
149 
151 {
152  if (cells_.empty())
153  {
154  createCells();
155  }
156 
157  return cells_;
158 }
159 
160 
162 {
163  if (patches_.empty())
164  {
165  createPatches();
166  }
167 
168  return patches_;
169 }
170 
171 
173 {
174  return topology().boundaryMesh().names();
175 }
176 
177 
178 //Foam::wordList Foam::blockMesh::patchTypes() const
179 //{
180 // return topology().boundaryMesh().types();
181 //}
182 //
183 //
184 //Foam::wordList Foam::blockMesh::patchPhysicalTypes() const
185 //{
186 // return topology().boundaryMesh().physicalTypes();
187 //}
188 
189 
191 {
192  label num = 0;
193 
194  forAll(*this, blocki)
195  {
196  if (operator[](blocki).zoneName().size())
197  {
198  num++;
199  }
200  }
201 
202  return num;
203 }
204 
205 
207 {
208  const pointField& pts = topology().points();
209 
210  forAll(pts, pI)
211  {
212  const point& pt = pts[pI];
213 
214  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
215  }
216 
217  const edgeList& edges = topology().edges();
218 
219  forAll(edges, eI)
220  {
221  const edge& e = edges[eI];
222 
223  os << "l " << e.start() + 1 << ' ' << e.end() + 1 << endl;
224  }
225 }
226 
227 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
wordList patchNames() const
Return patch names.
Definition: blockMesh.C:172
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:313
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
Foam::word regionName
defineDebugSwitch(blockMesh, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
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
A simple wrapper around bool so that it can be read as a word: true/false, on/off, yes/no, y/n, t/f, or none.
Definition: Switch.H:60
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:52
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1011
scalar scaleFactor() const
The scaling factor used to convert to metres.
Definition: blockMesh.C:133
pointField vertices(const blockVertexList &bvl)
stressControl lookup("compactNormalStress") >> compactNormalStress
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
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
wordList names() const
Return a list of patch names.
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
const Cmpt & x() const
Definition: VectorI.H:75
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
Define a block vertex.
Definition: blockVertex.H:48
~blockMesh()
Destructor.
Definition: blockMesh.C:83
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const pointField & vertices() const
Reference to point field defining the blockMesh.
Definition: blockMesh.C:97
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:150
Constant dispersed-phase particle diameter model.
label end() const
Return end vertex label.
Definition: edgeI.H:92
Input from memory buffer stream.
Definition: IStringStream.H:49
string str() const
Return the string.
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Output to memory buffer stream.
Definition: OStringStream.H:49
bool found
label start() const
Return start vertex label.
Definition: edgeI.H:81
Namespace for OpenFOAM.
label numZonedBlocks() const
Number of blocks with specified zones.
Definition: blockMesh.C:190