datToFoam.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-2016 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  datToFoam
29 
30 Group
31  grpMeshConversionUtilities
32 
33 Description
34  Reads in a datToFoam mesh file and outputs a points file.
35  Used in conjunction with blockMesh.
36 
37 \*---------------------------------------------------------------------------*/
38 
39 #include "argList.H"
40 #include "Time.H"
41 #include "Fstream.H"
42 #include "polyMesh.H"
43 #include "unitConversion.H"
44 
45 using namespace Foam;
46 
47 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
48 
49 int main(int argc, char *argv[])
50 {
52  (
53  "Reads in a datToFoam mesh file and outputs a points file.\n"
54  "Used in conjunction with blockMesh."
55  );
56 
58  argList::addArgument("dat file");
59 
60  argList args(argc, argv);
61 
62  if (!args.check())
63  {
64  FatalError.exit();
65  }
66 
67  #include "createTime.H"
68 
69  std::ifstream plot3dFile(args.get<fileName>(1));
70 
71  string line;
72  std::getline(plot3dFile, line);
73  std::getline(plot3dFile, line);
74 
75  IStringStream Istring(line);
76  word block;
77  string zoneName;
78  token punctuation;
79  label iPoints, jPoints;
80 
81  Istring >> block;
82  Istring >> block;
83  Istring >> zoneName;
84  Istring >> punctuation;
85  Istring >> block;
86  Istring >> iPoints;
87  Istring >> block;
88  Istring >> jPoints;
89 
90  Info<< "Number of vertices in i direction = " << iPoints << nl
91  << "Number of vertices in j direction = " << jPoints << endl;
92 
93  // We ignore the first layer of points in i and j the biconic meshes
94  const label nPointsij = (iPoints - 1)*(jPoints - 1);
95 
96  pointField points(nPointsij, Zero);
97 
98  for (direction comp = 0; comp < 2; ++comp)
99  {
100  label p = 0;
101 
102  for (label j = 0; j < jPoints; ++j)
103  {
104  for (label i = 0; i < iPoints; ++i)
105  {
106  scalar coord;
107  plot3dFile >> coord;
108 
109  // if statement ignores the first layer in i and j
110  if (i>0 && j>0)
111  {
112  points[p++][comp] = coord;
113  }
114  }
115  }
116  }
117 
118  // correct error in biconic meshes
119  for (point& p : points)
120  {
121  if (p.y() < 1e-07)
122  {
123  p.y() = 0;
124  }
125  }
126 
127  pointField pointsWedge(nPointsij*2, Zero);
128 
129  const scalar a = degToRad(0.1);
130  const tensor rotateZ
131  (
132  1.0, 0.0, 0.0,
133  0.0, 1.0, 0.0,
134  0.0, -::sin(a), ::cos(a)
135  );
136 
137  forAll(points, i)
138  {
139  pointsWedge[i] = (rotateZ & points[i]);
140  pointsWedge[i + nPointsij] =
142  (
143  vector(1.0, 1.0, -1.0),
144  pointsWedge[i]
145  );
146  }
147 
148  // Write (overwrite) points
149  {
150  IOobject iopoints
151  (
152  "points",
153  runTime.constant(), // instance
154  polyMesh::meshSubDir, // local
155  runTime
156  );
157 
158  if (!exists(iopoints.path()))
159  {
160  mkDir(iopoints.path());
161  }
162 
163  OFstream os(iopoints.objectPath(), runTime.writeStreamOption());
164 
165  Info<< "Writing points to: " << nl
166  << " " << os.name() << endl;
167 
168  iopoints.writeHeader(os, vectorIOField::typeName);
169  os << pointsWedge;
171  }
172 
173  Info<< "End\n" << endl;
174 
175  return 0;
176 }
177 
178 
179 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::Tensor< scalar >
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::block
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:58
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::exists
bool exists(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: MSwindows.C:625
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::cmptMultiply
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::OFstream::name
virtual const fileName & name() const
Read/write access to the name of the stream.
Definition: OSstream.H:107
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::IOobject::writeEndDivider
static Ostream & writeEndDivider(Ostream &os)
Write the standard end file divider.
Definition: IOobjectWriteHeader.C:142
Foam::sin
dimensionedScalar sin(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:264
Foam::polyMesh::meshSubDir
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:321
Foam::argList::addNote
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:412
unitConversion.H
Unit conversion functions.
Foam::argList
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:123
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::token
A token holds an item read from Istream.
Definition: token.H:68
polyMesh.H
Foam::argList::get
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:278
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::argList::addArgument
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
Definition: argList.C:301
Foam::Time::writeStreamOption
IOstreamOption writeStreamOption() const
The write stream option (format, compression, version)
Definition: Time.H:381
Foam::Field< vector >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
argList.H
Foam::FatalError
error FatalError
os
OBJstream os(runTime.globalPath()/outputName)
Foam::IStringStream
Input from string buffer, using a ISstream. Always UNCOMPRESSED.
Definition: StringStream.H:108
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::error::exit
void exit(const int errNo=1)
Exit : can be called for any error to exit program.
Definition: error.C:331
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::degToRad
constexpr scalar degToRad(const scalar deg) noexcept
Conversion from degrees to radians.
Definition: unitConversion.H:48
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Time.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Fstream.H
Foam::Vector< scalar >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::direction
uint8_t direction
Definition: direction.H:52
createTime.H
Foam::line
A line primitive.
Definition: line.H:53
Foam::argList::check
bool check(bool checkArgs=argList::argsMandatory(), bool checkOpts=true) const
Definition: argList.C:1908
Foam::argList::noParallel
static void noParallel()
Remove the parallel options.
Definition: argList.C:510
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Make a directory and return an error if it could not be created.
Definition: MSwindows.C:507
Foam::cos
dimensionedScalar cos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:265