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-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 \*---------------------------------------------------------------------------*/
25 
26 #include "blockMesh.H"
27 #include "Switch.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 bool Foam::blockMesh::blockMesh::verboseOutput(false);
32 
33 namespace Foam
34 {
35  defineDebugSwitch(blockMesh, 0);
36 }
37 
38 
39 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
40 
42 :
43  blockPointField_(dict.lookup("vertices")),
44  scaleFactor_(1.0),
45  topologyPtr_(createTopology(dict, regionName))
46 {
47  Switch fastMerge(dict.lookupOrDefault<Switch>("fastMerge", false));
48 
49  if (fastMerge)
50  {
52  }
53  else
54  {
55  calcMergeInfo();
56  }
57 }
58 
59 
60 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
61 
63 {
64  delete topologyPtr_;
65 }
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 void Foam::blockMesh::verbose(const bool on)
71 {
72  verboseOutput = on;
73 }
74 
75 
77 {
78  return blockPointField_;
79 }
80 
81 
83 {
84  if (!topologyPtr_)
85  {
87  << "topologyPtr_ not allocated"
88  << exit(FatalError);
89  }
90 
91  return *topologyPtr_;
92 }
93 
94 
96 {
97  const polyPatchList& patchTopologies = topology().boundaryMesh();
98 
99  PtrList<dictionary> patchDicts(patchTopologies.size());
100 
101  forAll(patchTopologies, patchI)
102  {
103  OStringStream os;
104  patchTopologies[patchI].write(os);
105  IStringStream is(os.str());
106  patchDicts.set(patchI, new dictionary(is));
107  }
108  return patchDicts;
109 }
110 
111 
112 Foam::scalar Foam::blockMesh::scaleFactor() const
113 {
114  return scaleFactor_;
115 }
116 
117 
119 {
120  if (points_.empty())
121  {
122  createPoints();
123  }
124 
125  return points_;
126 }
127 
128 
130 {
131  if (cells_.empty())
132  {
133  createCells();
134  }
135 
136  return cells_;
137 }
138 
139 
141 {
142  if (patches_.empty())
143  {
144  createPatches();
145  }
146 
147  return patches_;
148 }
149 
150 
152 {
153  return topology().boundaryMesh().names();
154 }
155 
156 
157 //Foam::wordList Foam::blockMesh::patchTypes() const
158 //{
159 // return topology().boundaryMesh().types();
160 //}
161 //
162 //
163 //Foam::wordList Foam::blockMesh::patchPhysicalTypes() const
164 //{
165 // return topology().boundaryMesh().physicalTypes();
166 //}
167 
168 
170 {
171  label num = 0;
172 
173  forAll(*this, blockI)
174  {
175  if (operator[](blockI).zoneName().size())
176  {
177  num++;
178  }
179  }
180 
181  return num;
182 }
183 
184 
186 {
187  const pointField& pts = topology().points();
188 
189  forAll(pts, pI)
190  {
191  const point& pt = pts[pI];
192 
193  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
194  }
195 
196  const edgeList& edges = topology().edges();
197 
198  forAll(edges, eI)
199  {
200  const edge& e = edges[eI];
201 
202  os << "l " << e.start() + 1 << ' ' << e.end() + 1 << endl;
203  }
204 }
205 
206 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
Foam::blockMesh::patchDicts
PtrList< dictionary > patchDicts() const
Get patch information from the topology mesh.
Definition: blockMesh.C:95
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:60
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
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
blockMesh.H
Foam::blockMesh::topology
const polyMesh & topology() const
Return the blockMesh topology as a polyMesh.
Definition: blockMesh.C:82
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
Foam::OStringStream::str
string str() const
Return the string.
Definition: OStringStream.H:107
patchDicts
PtrList< dictionary > patchDicts
Definition: readKivaGrid.H:537
Foam::blockMesh::writeTopology
void writeTopology(Ostream &) const
Writes edges of blockMesh in OBJ format.
Definition: blockMesh.C:185
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
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::defineDebugSwitch
defineDebugSwitch(blockMesh, 0)
Switch.H
Foam::blockMesh::~blockMesh
~blockMesh()
Destructor.
Definition: blockMesh.C:62
Foam::PtrList< Foam::dictionary >
Foam::blockMesh::scaleFactor
scalar scaleFactor() const
The scaling factor used to convert to metres.
Definition: blockMesh.C:112
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Vector::x
const Cmpt & x() const
Definition: VectorI.H:65
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::Vector::z
const Cmpt & z() const
Definition: VectorI.H:77
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
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::Vector< scalar >
Foam::List< cellShape >
Foam::blockMesh::patchNames
wordList patchNames() const
Return patch names.
Definition: blockMesh.C:151
Foam::OStringStream
Output to memory buffer stream.
Definition: OStringStream.H:49
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::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::Vector::y
const Cmpt & y() const
Definition: VectorI.H:71
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::blockMesh::cells
const cellShapeList & cells() const
Return cell shapes list.
Definition: blockMesh.C:129
Foam::blockMesh::calcMergeInfo
void calcMergeInfo()
Determine the merge info and the final number of cells/points.
Definition: blockMeshMerge.C:30
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress