All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
writeSTL.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-2019 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 #include "STLtriangle.H"
28 #include "primitivePatch.H"
29 #include "HashTable.H"
30 #include "hashSignedLabel.H"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 void Foam::triSurface::writeSTLASCII(const bool writeSorted, Ostream& os) const
35 {
37 
38  surfacePatchList patches(calcPatches(faceMap));
39 
40  if (writeSorted)
41  {
42  label faceIndex = 0;
44  {
45  // Print all faces belonging to this region
46  const surfacePatch& patch = patches[patchi];
47 
48  os << "solid " << patch.name() << endl;
49 
50  for
51  (
52  label patchFacei = 0;
53  patchFacei < patch.size();
54  patchFacei++
55  )
56  {
57  const label facei = faceMap[faceIndex++];
58 
59  const vector& n = faceNormals()[facei];
60 
61  os << " facet normal "
62  << n.x() << ' ' << n.y() << ' ' << n.z() << nl
63  << " outer loop" << endl;
64 
65  const labelledTri& f = (*this)[facei];
66  const point& pa = points()[f[0]];
67  const point& pb = points()[f[1]];
68  const point& pc = points()[f[2]];
69 
70  os << " vertex "
71  << pa.x() << ' ' << pa.y() << ' ' << pa.z() << nl
72  << " vertex "
73  << pb.x() << ' ' << pb.y() << ' ' << pb.z() << nl
74  << " vertex "
75  << pc.x() << ' ' << pc.y() << ' ' << pc.z() << nl
76  << " endloop" << nl
77  << " endfacet" << endl;
78  }
79 
80  os << "endsolid " << patch.name() << endl;
81  }
82  }
83  else
84  {
85  // Get patch (=compact region) per face
86  labelList patchIDs(size());
88  {
89  label facei = patches[patchi].start();
90 
91  forAll(patches[patchi], i)
92  {
93  patchIDs[faceMap[facei++]] = patchi;
94  }
95  }
96 
97  label currentPatchi = -1;
98 
99  forAll(*this, facei)
100  {
101  if (currentPatchi != patchIDs[facei])
102  {
103  if (currentPatchi != -1)
104  {
105  // Have already valid patch. Close it.
106  os << "endsolid " << patches[currentPatchi].name()
107  << nl;
108  }
109  currentPatchi = patchIDs[facei];
110  os << "solid " << patches[currentPatchi].name() << nl;
111  }
112 
113  const vector& n = faceNormals()[facei];
114 
115  os << " facet normal "
116  << n.x() << ' ' << n.y() << ' ' << n.z() << nl
117  << " outer loop" << endl;
118 
119  const labelledTri& f = (*this)[facei];
120  const point& pa = points()[f[0]];
121  const point& pb = points()[f[1]];
122  const point& pc = points()[f[2]];
123 
124  os << " vertex "
125  << pa.x() << ' ' << pa.y() << ' ' << pa.z() << nl
126  << " vertex "
127  << pb.x() << ' ' << pb.y() << ' ' << pb.z() << nl
128  << " vertex "
129  << pc.x() << ' ' << pc.y() << ' ' << pc.z() << nl
130  << " endloop" << nl
131  << " endfacet" << endl;
132  }
133 
134  if (currentPatchi != -1)
135  {
136  os << "endsolid " << patches[currentPatchi].name()
137  << nl;
138  }
139  }
140 }
141 
142 
143 void Foam::triSurface::writeSTLBINARY(std::ostream& os) const
144 {
145  // Write the STL header
146  string header("Foam binary STL");
147  header.resize(STLheaderSize);
148  os.write(header.c_str(), STLheaderSize);
149 
150  label nTris = size();
151  os.write(reinterpret_cast<char*>(&nTris), sizeof(unsigned int));
152 
153  const vectorField& normals = faceNormals();
154 
155  forAll(*this, facei)
156  {
157  const labelledTri& f = (*this)[facei];
158 
159  // Convert vector into STL single precision
160  STLpoint n(normals[facei]);
161  STLpoint pa(points()[f[0]]);
162  STLpoint pb(points()[f[1]]);
163  STLpoint pc(points()[f[2]]);
164 
165  STLtriangle stlTri(n, pa, pb, pc, f.region());
166 
167  stlTri.write(os);
168  }
169 }
170 
171 
172 // ************************************************************************* //
#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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
List< surfacePatch > surfacePatchList
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:294
const Field< PointType > & points() const
Return reference to global points.
List< label > labelList
A List of labels.
Definition: labelList.H:56
const Field< PointType > & faceNormals() const
Return face normals for patch.
static const char nl
Definition: Ostream.H:260
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
label patchi
vector point
Point is a vector.
Definition: point.H:41
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
label size() const
Return the number of elements in the UList.
Definition: ListI.H:171