writeOFF.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::writeOFF(const bool writeSorted, Ostream& os) const
36 {
37  // Write header
38  os << "OFF" << endl
39  << "# Geomview OFF file" << endl
40  << "# Regions:" << endl;
41 
43  surfacePatchList patches(calcPatches(faceMap));
44 
45  // Print patch names as comment
47  {
48  os << "# " << patchi << " "
49  << patches[patchi].name() << endl;
50  }
51  os << nl << endl;
52 
53  const pointField& ps = points();
54 
55  os << "# nPoints nTriangles nEdges" << endl
56  << ps.size()
57  << ' ' << size()
58  << ' ' << nEdges()
59  << nl << endl;
60 
61  // Write vertex coords
62  forAll(ps, pointi)
63  {
64  os << ps[pointi].x() << ' '
65  << ps[pointi].y() << ' '
66  << ps[pointi].z() << " #" << pointi << endl;
67  }
68 
69  os << endl;
70 
71  if (writeSorted)
72  {
73  label faceIndex = 0;
74 
76  {
77  // Print all faces belonging to this patch
78 
79  for
80  (
81  label patchFacei = 0;
82  patchFacei < patches[patchi].size();
83  patchFacei++
84  )
85  {
86  const label facei = faceMap[faceIndex++];
87 
88  os << "3 "
89  << operator[](facei)[0] << ' '
90  << operator[](facei)[1] << ' '
91  << operator[](facei)[2] << ' '
92  << operator[](facei).region()
93  << endl;
94  }
95  }
96  }
97  else
98  {
99  forAll(*this, facei)
100  {
101  os << "3 "
102  << operator[](facei)[0] << ' '
103  << operator[](facei)[1] << ' '
104  << operator[](facei)[2] << ' '
105  << operator[](facei).region()
106  << endl;
107  }
108  }
109 }
110 
111 
112 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
113 
114 } // End namespace Foam
115 
116 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nEdges() const
Return number of edges in patch.
const Field< PointType > & points() const
Return reference to global points.
T & operator[](const label)
Return element of UList.
Definition: UListI.H:167
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
label patchi
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
static const char nl
Definition: Ostream.H:266
List< surfacePatch > surfacePatchList