All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
writeGTS.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-2018 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 "triSurface.H"
27 
28 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
34 
35 void triSurface::writeGTS(const bool writeSorted, Ostream& os) const
36 {
37  // Write header
38  os << "# GTS file" << endl
39  << "# Regions:" << endl;
40 
42 
43  surfacePatchList patches(calcPatches(faceMap));
44 
45  // Print patch names as comment
47  {
48  os << "# " << patchi << " "
49  << patches[patchi].name() << endl;
50  }
51  os << "#" << endl;
52 
53 
54  const pointField& ps = points();
55 
56  os << "# nPoints nEdges nTriangles" << endl
57  << ps.size() << ' ' << nEdges() << ' ' << size() << endl;
58 
59  // Write vertex coords
60 
61  forAll(ps, pointi)
62  {
63  os << ps[pointi].x() << ' '
64  << ps[pointi].y() << ' '
65  << ps[pointi].z() << endl;
66  }
67 
68  // Write edges.
69  // Note: edges are in local point labels so convert
70  const edgeList& es = edges();
71  const labelList& meshPts = meshPoints();
72 
73  forAll(es, edgei)
74  {
75  os << meshPts[es[edgei].start()] + 1 << ' '
76  << meshPts[es[edgei].end()] + 1 << endl;
77  }
78 
79  // Write faces in terms of edges.
80  const labelListList& faceEs = faceEdges();
81 
82  if (writeSorted)
83  {
84  label faceIndex = 0;
86  {
87  for
88  (
89  label patchFacei = 0;
90  patchFacei < patches[patchi].size();
91  patchFacei++
92  )
93  {
94  const label facei = faceMap[faceIndex++];
95 
96  const labelList& fEdges = faceEdges()[facei];
97 
98  os << fEdges[0] + 1 << ' '
99  << fEdges[1] + 1 << ' '
100  << fEdges[2] + 1 << ' '
101  << (*this)[facei].region() << endl;
102  }
103  }
104  }
105  else
106  {
107  forAll(faceEs, facei)
108  {
109  const labelList& fEdges = faceEdges()[facei];
110 
111  os << fEdges[0] + 1 << ' '
112  << fEdges[1] + 1 << ' '
113  << fEdges[2] + 1 << ' '
114  << (*this)[facei].region() << endl;
115  }
116  }
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
121 
122 } // End namespace Foam
123 
124 // ************************************************************************* //
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
List< surfacePatch > surfacePatchList
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
List< edge > edgeList
Definition: edgeList.H:38
const Field< PointType > & points() const
Return reference to global points.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< label > labelList
A List of labels.
Definition: labelList.H:56
label nEdges() const
Return number of edges in patch.
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
label patchi
const labelListList & faceEdges() const
Return face-edge addressing.
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
Namespace for OpenFOAM.