triSurfaceGTSformat.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 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 \*---------------------------------------------------------------------------*/
28 
29 #include "triSurface.H"
30 #include "OFstream.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 void Foam::triSurface::writeGTS
35 (
36  const fileName& filename,
37  const bool sort
38 ) const
39 {
40  OFstream os(filename);
41  if (!os.good())
42  {
44  << "Cannot write file " << filename << nl
45  << exit(FatalError);
46  }
47 
48  // Write header
49  os << "# GTS file" << endl
50  << "# Regions:" << endl;
51 
53  surfacePatchList patches(calcPatches(faceMap));
54 
55  // Print patch names as comment
56  forAll(patches, patchi)
57  {
58  os << "# " << patchi << " "
59  << patches[patchi].name() << nl;
60  }
61  os << "#" << nl;
62 
63  const pointField& pts = points();
64 
65  os << "# nPoints nEdges nTriangles" << nl
66  << pts.size() << ' ' << nEdges() << ' ' << size() << nl;
67 
68  // Write vertex coords
69  for (const point& pt : pts)
70  {
71  os << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
72  }
73 
74  // Write edges.
75  // Note: edges are in local point labels so convert
76  const edgeList& es = edges();
77  const labelList& meshPts = meshPoints();
78 
79  for (const edge& e : es)
80  {
81  os << meshPts[e.start()] + 1 << ' '
82  << meshPts[e.end()] + 1 << nl;
83  }
84 
85  // Write faces in terms of edges.
86  const labelListList& faceEs = faceEdges();
87 
88  if (sort)
89  {
90  label faceIndex = 0;
91  for (const surfacePatch& p : patches)
92  {
93  for (label nLocal = p.size(); nLocal--; ++faceIndex)
94  {
95  const label facei = faceMap[faceIndex];
96 
97  const labelList& fEdges = faceEdges()[facei];
98 
99  os << fEdges[0] + 1 << ' '
100  << fEdges[1] + 1 << ' '
101  << fEdges[2] + 1 << ' '
102  << (*this)[facei].region() << nl;
103  }
104  }
105  }
106  else
107  {
108  forAll(faceEs, facei)
109  {
110  const labelList& fEdges = faceEdges()[facei];
111 
112  os << fEdges[0] + 1 << ' '
113  << fEdges[1] + 1 << ' '
114  << fEdges[2] + 1 << ' '
115  << (*this)[facei].region() << nl;
116  }
117  }
118 }
119 
120 
121 // ************************************************************************* //
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const geometricSurfacePatchList & patches() const noexcept
Definition: triSurface.H:509
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Field< point_type > & points() const noexcept
Return reference to global points.
label size() const noexcept
The number of elements in the container.
Definition: UList.H:671
OBJstream os(runTime.globalPath()/outputName)
vector point
Point is a vector.
Definition: point.H:37
const labelListList & faceEdges() const
Return face-edge addressing.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
List< surfacePatch > surfacePatchList
List of surfacePatch.
const pointField & pts