plot3dToFoam.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  plot3dToFoam
29 
30 Group
31  grpMeshConversionUtilities
32 
33 Description
34  Plot3d mesh (ascii/formatted format) converter.
35 
36  Work in progress! Handles ascii multiblock (and optionally singleBlock)
37  format.
38  By default expects blanking. Use -noBlank if none.
39  Use -2D \a thickness if 2D.
40 
41  Niklas Nordin has experienced a problem with lefthandedness of the blocks.
42  The code should detect this automatically - see hexBlock::readPoints but
43  if this goes wrong just set the blockHandedness_ variable to 'right'
44  always.
45 
46 \*---------------------------------------------------------------------------*/
47 
48 #include "argList.H"
49 #include "Time.H"
50 #include "IFstream.H"
51 #include "hexBlock.H"
52 #include "polyMesh.H"
53 #include "wallPolyPatch.H"
54 #include "symmetryPolyPatch.H"
55 #include "cellShape.H"
56 #include "mergePoints.H"
57 
58 using namespace Foam;
59 
60 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
61 
62 int main(int argc, char *argv[])
63 {
65  (
66  "Plot3d mesh (ascii/formatted format) converter"
67  );
69  argList::addArgument("PLOT3D geom file");
71  (
72  "scale",
73  "factor",
74  "Geometry scaling factor - default is 1"
75  );
77  (
78  "noBlank",
79  "Skip blank items"
80  );
82  (
83  "singleBlock",
84  "Input is a single block"
85  );
87  (
88  "2D",
89  "thickness",
90  "Use when converting a 2-D mesh (applied before scale)"
91  );
92 
93  argList args(argc, argv);
94 
95  if (!args.check())
96  {
97  FatalError.exit();
98  }
99 
100  const scalar scaleFactor = args.getOrDefault<scalar>("scale", 1);
101 
102  const bool readBlank = !args.found("noBlank");
103  const bool singleBlock = args.found("singleBlock");
104  scalar twoDThickness = -1;
105  if (args.readIfPresent("2D", twoDThickness))
106  {
107  Info<< "Reading 2D case by extruding points by " << twoDThickness
108  << " in z direction." << nl << endl;
109  }
110 
111 
112  #include "createTime.H"
113 
114  IFstream plot3dFile(args.get<fileName>(1));
115 
116  // Read the plot3d information using a fixed format reader.
117  // Comments in the file are in C++ style, so the stream parser will remove
118  // them with no intervention
119  label nblock;
120 
121  if (singleBlock)
122  {
123  nblock = 1;
124  }
125  else
126  {
127  plot3dFile >> nblock;
128  }
129 
130  Info<< "Reading " << nblock << " blocks" << endl;
131 
132  PtrList<hexBlock> blocks(nblock);
133 
134  {
135  label nx, ny, nz;
136 
137  forAll(blocks, blockI)
138  {
139  if (twoDThickness > 0)
140  {
141  // Fake second set of points (done in readPoints below)
142  plot3dFile >> nx >> ny;
143  nz = 2;
144  }
145  else
146  {
147  plot3dFile >> nx >> ny >> nz;
148  }
149 
150  Info<< "block " << blockI << " nx:" << nx
151  << " ny:" << ny << " nz:" << nz << endl;
152 
153  blocks.set(blockI, new hexBlock(nx, ny, nz));
154  }
155  }
156 
157  Info<< "Reading block points" << endl;
158  label sumPoints(0);
159  label nMeshCells(0);
160 
161  forAll(blocks, blockI)
162  {
163  Info<< "block " << blockI << ":" << nl;
164  blocks[blockI].readPoints(readBlank, twoDThickness, plot3dFile);
165  sumPoints += blocks[blockI].nBlockPoints();
166  nMeshCells += blocks[blockI].nBlockCells();
167  Info<< nl;
168  }
169 
170  pointField points(sumPoints);
171  labelList blockOffsets(blocks.size());
172  sumPoints = 0;
173  forAll(blocks, blockI)
174  {
175  const pointField& blockPoints = blocks[blockI].points();
176  blockOffsets[blockI] = sumPoints;
177  forAll(blockPoints, i)
178  {
179  points[sumPoints++] = blockPoints[i];
180  }
181  }
182 
183  // From old to new master point
184  labelList oldToNew;
185  pointField newPoints;
186 
187  // Merge points
189  (
190  points,
191  SMALL,
192  false,
193  oldToNew,
194  newPoints
195  );
196 
197  Info<< "Merged points within " << SMALL << " distance. Merged from "
198  << oldToNew.size() << " down to " << newPoints.size()
199  << " points." << endl;
200 
201  // Scale the points
202  if (scaleFactor > 1.0 + SMALL || scaleFactor < 1.0 - SMALL)
203  {
204  newPoints *= scaleFactor;
205  }
206 
207  Info<< "Creating cells" << endl;
208 
209  cellShapeList cellShapes(nMeshCells);
210 
212 
213  label nCreatedCells = 0;
214 
215  forAll(blocks, blockI)
216  {
217  labelListList curBlockCells = blocks[blockI].blockCells();
218 
219  forAll(curBlockCells, blockCelli)
220  {
221  labelList cellPoints(curBlockCells[blockCelli].size());
222 
223  forAll(cellPoints, pointi)
224  {
225  cellPoints[pointi] =
226  oldToNew
227  [
228  curBlockCells[blockCelli][pointi]
229  + blockOffsets[blockI]
230  ];
231  }
232 
233  // Do automatic collapse from hex.
234  cellShapes[nCreatedCells].reset(hex, cellPoints, true);
235 
236  nCreatedCells++;
237  }
238  }
239 
240  Info<< "Creating boundary patches" << endl;
241 
243  wordList patchNames(0);
244  wordList patchTypes(0);
245  word defaultFacesName = "defaultFaces";
246  word defaultFacesType = wallPolyPatch::typeName;
247  wordList patchPhysicalTypes(0);
248 
250  (
251  IOobject
252  (
254  runTime.constant(),
255  runTime
256  ),
257  std::move(newPoints),
258  cellShapes,
259  boundary,
260  patchNames,
261  patchTypes,
264  patchPhysicalTypes
265  );
266 
267  // Set the precision of the points data to 10
269 
270  Info<< "Writing polyMesh" << endl;
272  pShapeMesh.write();
273 
274  Info<< "End\n" << endl;
275 
276  return 0;
277 }
278 
279 
280 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::polyMesh::defaultRegion
static word defaultRegion
Definition: polyMesh.H:314
Foam::argList::getOrDefault
T getOrDefault(const word &optName, const T &deflt) const
Definition: argListI.H:300
Foam::IFstream
Input from file stream, using an ISstream.
Definition: IFstream.H:49
Foam::cellModel::HEX
@ HEX
hex
Definition: cellModel.H:77
wallPolyPatch.H
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
symmetryPolyPatch.H
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:119
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
polyMesh.H
defaultFacesType
word defaultFacesType
Definition: readKivaGrid.H:456
Foam::argList::get
T get(const label index) const
Definition: argListI.H:271
Foam::argList::readIfPresent
bool readIfPresent(const word &optName, T &val) const
Definition: argListI.H:316
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
patchTypes
wordList patchTypes(nPatches)
Foam::argList::addArgument
static void addArgument(const string &argName, const string &usage="")
Definition: argList.C:294
Foam::Field
Generic templated field type.
Definition: Field.H:59
Foam::regIOobject::write
virtual bool write(const bool valid=true) const
Definition: regIOobjectWrite.C:125
Foam::Info
messageStream Info
pShapeMesh
polyMesh pShapeMesh(IOobject(polyMesh::defaultRegion, runTime.constant(), runTime), std::move(points), cellShapes, boundary, patchNames, patchDicts, defaultFacesName, defaultFacesType)
Foam::polyMesh::removeFiles
void removeFiles(const fileName &instanceDir) const
Definition: polyMesh.C:1318
argList.H
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:55
patchNames
wordList patchNames(nPatches)
Foam::cellModel::ref
static const cellModel & ref(const modelType model)
Definition: cellModels.C:150
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Definition: hashSets.C:40
IFstream.H
Foam::FatalError
error FatalError
Foam
Definition: atmBoundaryLayer.C:26
Foam::error::exit
void exit(const int errNo=1)
Definition: error.C:324
Foam::mergePoints
label mergePoints(const PointList &points, const scalar mergeTol, const bool verbose, labelList &pointMap, typename PointList::const_reference origin=PointList::value_type::zero)
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:446
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Definition: argList.C:317
Foam::IOstream::defaultPrecision
static unsigned int defaultPrecision() noexcept
Definition: IOstream.H:338
Time.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
defaultFacesName
word defaultFacesName
Definition: readKivaGrid.H:455
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
points
const pointField & points
Definition: gmvOutputHeader.H:1
createTime.H
mergePoints.H
Merge points. See below.
cellShape.H
cellShapes
cellShapeList cellShapes
Definition: createBlockMesh.H:3
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:68
Foam::hexBlock
Hex block definition used in the cfx converter.
Definition: hexBlock.H:48
Foam::argList::check
bool check(bool checkArgs=argList::argsMandatory(), bool checkOpts=true) const
Definition: argList.C:1901
Foam::argList::noParallel
static void noParallel()
Definition: argList.C:503
Foam::TimePaths::constant
const word & constant() const
Definition: TimePathsI.H:89
Foam::argList::addOption
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Definition: argList.C:328
args
Foam::argList args(argc, argv)
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171