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 | 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 Application
25  datToFoam
26 
27 Description
28  Reads in a datToFoam mesh file and outputs a points file. Used in
29  conjunction with blockMesh.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "argList.H"
34 #include "Time.H"
35 #include "IFstream.H"
36 #include "OFstream.H"
37 #include "pointField.H"
38 #include "unitConversion.H"
39 
40 using namespace Foam;
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
47  argList::validArgs.append("dat file");
48 
49  argList args(argc, argv);
50 
51  if (!args.check())
52  {
53  FatalError.exit();
54  }
55 
56  #include "createTime.H"
57 
58  std::ifstream plot3dFile(args.args()[1].c_str());
59 
60  string line;
61  std::getline(plot3dFile, line);
62  std::getline(plot3dFile, line);
63 
64  IStringStream Istring(line);
65  word block;
66  string zoneName;
67  token punctuation;
68  label iPoints;
69  label jPoints;
70 
71  Istring >> block;
72  Istring >> block;
73  Istring >> zoneName;
74  Istring >> punctuation;
75  Istring >> block;
76  Istring >> iPoints;
77  Istring >> block;
78  Istring >> jPoints;
79 
80  Info<< "Number of vertices in i direction = " << iPoints << endl
81  << "Number of vertices in j direction = " << jPoints << endl;
82 
83  // We ignore the first layer of points in i and j the biconic meshes
84  label nPointsij = (iPoints - 1)*(jPoints - 1);
85 
86  pointField points(nPointsij, Zero);
87 
88  for (direction comp = 0; comp < 2; comp++)
89  {
90  label p(0);
91 
92  for (label j = 0; j < jPoints; j++)
93  {
94  for (label i = 0; i < iPoints; i++)
95  {
96  double coord;
97  plot3dFile >> coord;
98 
99  // if statement ignores the first layer in i and j
100  if (i>0 && j>0)
101  {
102  points[p++][comp] = coord;
103  }
104  }
105  }
106  }
107 
108  // correct error in biconic meshes
109  forAll(points, i)
110  {
111  if (points[i][1] < 1e-07)
112  {
113  points[i][1] = 0.0;
114  }
115  }
116 
117  pointField pointsWedge(nPointsij*2, Zero);
118 
119  fileName pointsFile(runTime.constantPath()/"points.tmp");
120  OFstream pFile(pointsFile);
121 
122  scalar a(degToRad(0.1));
123  tensor rotateZ =
124  tensor
125  (
126  1.0, 0.0, 0.0,
127  0.0, 1.0, 0.0,
128  0.0, -::sin(a), ::cos(a)
129  );
130 
131  forAll(points, i)
132  {
133  pointsWedge[i] = (rotateZ & points[i]);
134  pointsWedge[i+nPointsij] = cmptMultiply
135  (
136  vector(1.0, 1.0, -1.0),
137  pointsWedge[i]
138  );
139  }
140 
141  Info<< "Writing points to: " << nl
142  << " " << pointsFile << endl;
143  pFile << pointsWedge;
144 
145  Info<< "End" << endl;
146 
147  return 0;
148 }
149 
150 
151 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
A class for handling file names.
Definition: fileName.H:69
const double e
Elementary charge.
Definition: doubleFloat.H:78
error FatalError
Unit conversion functions.
Output to file stream.
Definition: OFstream.H:82
uint8_t direction
Definition: direction.H:45
A token holds items read from Istream.
Definition: token.H:69
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void noParallel()
Remove the parallel options.
Definition: argList.C:148
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:154
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void exit(const int errNo=1)
Exit : can be called for any error to exit program.
Definition: error.C:168
scalar degToRad(const scalar deg)
Conversion from degrees to radians.
dimensionedScalar cos(const dimensionedScalar &ds)
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
Extract command arguments and options from the supplied argc and argv parameters. ...
Definition: argList.H:102
static const zero Zero
Definition: zero.H:91
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar sin(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:262
Input from memory buffer stream.
Definition: IStringStream.H:49
messageStream Info
bool check(bool checkArgs=true, bool checkOpts=true) const
Check argument list.
Definition: argList.C:1246
volScalarField & p
Foam::argList args(argc, argv)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Namespace for OpenFOAM.
const stringList & args() const
Return arguments.
Definition: argListI.H:66