datToFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 
39 using namespace Foam;
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
46  argList::validArgs.append("dat file");
47 
48  argList args(argc, argv);
49 
50  if (!args.check())
51  {
52  FatalError.exit();
53  }
54 
55  #include "createTime.H"
56 
57  std::ifstream plot3dFile(args.args()[1].c_str());
58 
59  string line;
60  std::getline(plot3dFile, line);
61  std::getline(plot3dFile, line);
62 
63  IStringStream Istring(line);
64  word block;
65  string zoneName;
66  token punctuation;
67  label iPoints;
68  label jPoints;
69 
70  Istring >> block;
71  Istring >> block;
72  Istring >> zoneName;
73  Istring >> punctuation;
74  Istring >> block;
75  Istring >> iPoints;
76  Istring >> block;
77  Istring >> jPoints;
78 
79  Info<< "Number of vertices in i direction = " << iPoints << endl
80  << "Number of vertices in j direction = " << jPoints << endl;
81 
82  // We ignore the first layer of points in i and j the biconic meshes
83  label nPointsij = (iPoints - 1)*(jPoints - 1);
84 
85  pointField points(nPointsij, Zero);
86 
87  for (direction comp = 0; comp < 2; comp++)
88  {
89  label p(0);
90 
91  for (label j = 0; j < jPoints; j++)
92  {
93  for (label i = 0; i < iPoints; i++)
94  {
95  double coord;
96  plot3dFile >> coord;
97 
98  // if statement ignores the first layer in i and j
99  if (i>0 && j>0)
100  {
101  points[p++][comp] = coord;
102  }
103  }
104  }
105  }
106 
107  // correct error in biconic meshes
108  forAll(points, i)
109  {
110  if (points[i][1] < 1e-07)
111  {
112  points[i][1] = 0.0;
113  }
114  }
115 
116  pointField pointsWedge(nPointsij*2, Zero);
117 
118  fileName pointsFile(runTime.path()/runTime.constant()/"points.tmp");
119  OFstream pFile(pointsFile);
120 
121  scalar a(degToRad(0.1));
122  tensor rotateZ =
123  tensor
124  (
125  1.0, 0.0, 0.0,
126  0.0, 1.0, 0.0,
127  0.0, -::sin(a), ::cos(a)
128  );
129 
130  forAll(points, i)
131  {
132  pointsWedge[i] = (rotateZ & points[i]);
133  pointsWedge[i+nPointsij] = cmptMultiply
134  (
135  vector(1.0, 1.0, -1.0),
136  pointsWedge[i]
137  );
138  }
139 
140  Info<< "Writing points to: " << nl
141  << " " << pointsFile << endl;
142  pFile << pointsWedge;
143 
144  Info<< "End" << endl;
145 
146  return 0;
147 }
148 
149 
150 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Input from memory buffer stream.
Definition: IStringStream.H:52
Output to file stream.
Definition: OFstream.H:86
Extract command arguments and options from the supplied argc and argv parameters.
Definition: argList.H:103
const stringList & args() const
Return arguments.
Definition: argListI.H:72
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static SLList< string > validArgs
A list of valid (mandatory) arguments.
Definition: argList.H:153
bool check(bool checkArgs=true, bool checkOpts=true) const
Check argument list.
Definition: argList.C:1400
Creates a single block of cells from point coordinates, numbers of cells in each direction and an exp...
Definition: block.H:66
void exit(const int errNo=1)
Exit : can be called for any error to exit program.
Definition: error.C:125
A class for handling file names.
Definition: fileName.H:82
A line primitive.
Definition: line.H:71
A token holds items read from Istream.
Definition: token.H:73
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
const pointField & points
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
scalar degToRad(const scalar deg)
Convert degrees to radians.
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
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
static const char nl
Definition: Ostream.H:266
dimensionedScalar cos(const dimensionedScalar &ds)
uint8_t direction
Definition: direction.H:45
Foam::argList args(argc, argv)
volScalarField & p