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 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "triSurface.H"
29 #include "OFstream.H"
30 
31 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
32 
33 void Foam::triSurface::writeGTS
34 (
35  const fileName& filename,
36  const bool sort
37 ) const
38 {
39  OFstream os(filename);
40  if (!os.good())
41  {
43  << "Cannot open file for writing " << filename
44  << exit(FatalError);
45  }
46 
47  // Write header
48  os << "# GTS file" << endl
49  << "# Regions:" << endl;
50 
52  surfacePatchList patches(calcPatches(faceMap));
53 
54  // Print patch names as comment
55  forAll(patches, patchi)
56  {
57  os << "# " << patchi << " "
58  << patches[patchi].name() << nl;
59  }
60  os << "#" << nl;
61 
62  const pointField& pts = points();
63 
64  os << "# nPoints nEdges nTriangles" << nl
65  << pts.size() << ' ' << nEdges() << ' ' << size() << nl;
66 
67  // Write vertex coords
68  for (const point& pt : pts)
69  {
70  os << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
71  }
72 
73  // Write edges.
74  // Note: edges are in local point labels so convert
75  const edgeList& es = edges();
76  const labelList& meshPts = meshPoints();
77 
78  for (const edge& e : es)
79  {
80  os << meshPts[e.start()] + 1 << ' '
81  << meshPts[e.end()] + 1 << nl;
82  }
83 
84  // Write faces in terms of edges.
85  const labelListList& faceEs = faceEdges();
86 
87  if (sort)
88  {
89  label faceIndex = 0;
90  for (const surfacePatch& p : patches)
91  {
92  const label nLocalFaces = p.size();
93 
94  for (label i = 0; i<nLocalFaces; ++i)
95  {
96  const label facei = faceMap[faceIndex++];
97 
98  const labelList& fEdges = faceEdges()[facei];
99 
100  os << fEdges[0] + 1 << ' '
101  << fEdges[1] + 1 << ' '
102  << fEdges[2] + 1 << ' '
103  << (*this)[facei].region() << nl;
104  }
105  }
106  }
107  else
108  {
109  forAll(faceEs, facei)
110  {
111  const labelList& fEdges = faceEdges()[facei];
112 
113  os << fEdges[0] + 1 << ' '
114  << fEdges[1] + 1 << ' '
115  << fEdges[2] + 1 << ' '
116  << (*this)[facei].region() << nl;
117  }
118  }
119 }
120 
121 
122 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:94
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::points
const Field< point > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
Foam::edgeList
List< edge > edgeList
A List of edges.
Definition: edgeList.H:63
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:238
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatch.H:317
Foam::surfacePatchList
List< surfacePatch > surfacePatchList
Definition: surfacePatchList.H:46
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::faceEdges
const labelListList & faceEdges() const
Return face-edge addressing.
Definition: PrimitivePatch.C:338
triSurface.H
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
OFstream.H
Foam::label
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:62
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:241
Foam::FatalError
error FatalError
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418