iRoCS Toolbox  1.1.0
STLFileWriter.hh
Go to the documentation of this file.
1 /**************************************************************************
2  *
3  * Copyright (C) 2015 Thorsten Falk
4  *
5  * Image Analysis Lab, University of Freiburg, Germany
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 3 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, write to the Free Software Foundation,
19  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
20  *
21  **************************************************************************/
22 
23 #ifndef STLFILEWRITER_HH
24 #define STLFILEWRITER_HH
25 
26 #ifdef HAVE_CONFIG_H
27 #include <config.hh>
28 #endif
29 
30 #include <fstream>
31 #include <vector>
32 #include <blitz/array.h>
33 
34 namespace atb
35 {
36 
38  {
39  std::ofstream _outFile;
40  std::string _objectName;
41  float _tailRadius;
42  float _headRadius;
43  int _nSegments;
44  float _borderRadius;
45  int _nEllipseSegments;
46 
47  public:
48 
49  STLFileWriter(const std::string& fileName,
50  const std::string& objectName="dummy")
51  : _outFile(fileName.c_str()),
52  _objectName(objectName),
53  _tailRadius(0.5),
54  _headRadius(2),
55  _nSegments(20),
56  _borderRadius(0.5),
57  _nEllipseSegments(72)
58  {
59  _outFile << "solid " << _objectName << "\n";
60  }
61 
63  {
64  _outFile << "endsolid " << _objectName << "\n";
65  }
66 
67 
68 
69  void drawTriangle(const blitz::TinyVector<float,3>& p1,
70  const blitz::TinyVector<float,3>& p2,
71  const blitz::TinyVector<float,3>& p3,
72  const blitz::TinyVector<float,3>& normal =
73  blitz::TinyVector<float,3>(0,0,0))
74  {
75  _outFile << "facet normal "
76  << normal(0) << " " << normal(1) << " " << normal(2)
77  << "\n"
78  << " outer loop\n"
79  << " vertex "
80  << p1(0) << " " << p1(1) << " " << p1(2) << "\n"
81  << " vertex "
82  << p2(0) << " " << p2(1) << " " << p2(2) << "\n"
83  << " vertex "
84  << p3(0) << " " << p3(1) << " " << p3(2) << "\n"
85  << " endloop\n"
86  << "endfacet\n";
87  }
88 
90  const blitz::TinyVector<float,3>& origin,
91  const blitz::TinyVector<float,3>& nx,
92  const blitz::TinyVector<float,3>& ny,
93  const blitz::TinyVector<float,3>& nz,
94  const std::vector< blitz::TinyVector<float,2> >& contour,
95  const std::vector< blitz::TinyVector<float,2> >& normals,
96  const int nAngles)
97  {
98  // make a rotation object from this contour
99  for( int i = 0; i < nAngles; ++i)
100  {
101  float phi1 = float(i)/nAngles*2*M_PI;
102  float phi2 = float(i+1)/nAngles*2*M_PI;
103  blitz::TinyVector<float,3> rotNy1 =
104  sin( phi1) * ny + cos( phi1) * nz;
105 
106  blitz::TinyVector<float,3> rotNy2 =
107  sin( phi2) * ny + cos( phi2) * nz;
108 
109  blitz::TinyVector<float,3> rotNyMean =
110  (rotNy1 + rotNy2)/2;
111  rotNyMean /= sqrt( blitz::dot( rotNyMean, rotNyMean));
112 
113  for( size_t j = 0; j < contour.size()-1; ++j)
114  {
115  blitz::TinyVector<float,3> p1 =
116  origin
117  + nx * contour[j](0)
118  + rotNy1 * contour[j](1);
119 
120  blitz::TinyVector<float,3> p2 =
121  origin
122  + nx * contour[j](0)
123  + rotNy2 * contour[j](1);
124 
125  blitz::TinyVector<float,3> p3 =
126  origin
127  + nx * contour[j+1](0)
128  + rotNy1 * contour[j+1](1);
129 
130  blitz::TinyVector<float,3> p4 =
131  origin
132  + nx * contour[j+1](0)
133  + rotNy2 * contour[j+1](1);
134 
135  blitz::TinyVector<float,3> normal =
136  nx * normals[j](0)
137  + rotNyMean * normals[j](1);
138 
139  // triangle 1 (p1 -- p2 -- p3)
140  if(blitz::any(p1 != p2))
141  {
142  drawTriangle( p1, p2, p3, normal);
143  }
144 
145  // triangle 2 (p2 -- p4 -- p3)
146  if(blitz::any(p3 != p4))
147  {
148  drawTriangle( p2, p4, p3, normal);
149  }
150  }
151 
152  }
153  }
154 
155 
156  void drawCylinder(const blitz::TinyVector<float,3>& origin,
157  const blitz::TinyVector<float,3>& dxyz,
158  double radius = -1)
159  {
160  if (radius <= 0) radius = _tailRadius;
161 
162  // compute x- and y- unit vectors
163  float length = sqrt( blitz::dot( dxyz, dxyz));
164 
165  blitz::TinyVector<float,3> nx = dxyz / length;
166 
167  blitz::TinyVector<float,3> unitCol(0,0,1);
168  blitz::TinyVector<float,3> ny;
169  blitz::TinyVector<float,3> nz;
170 
171  ny = blitz::cross( nx, unitCol);
172  nz = blitz::cross( nx, ny);
173 
174  ny /= sqrt( blitz::dot( ny, ny));
175  nz /= sqrt( blitz::dot( nz, nz));
176 
177  // create the contour of the arrow
178  std::vector<blitz::TinyVector<float,2> > contour(4);
179  contour[0] = length, 0;
180  contour[1] = length, radius;
181  contour[2] = 0, radius;
182  contour[3] = 0,0;
183 
184  std::vector<blitz::TinyVector<float,2> > normals(4);
185  normals[0] = 1, 0;
186  normals[1] = 0, 1;
187  normals[2] = 0, 1;
188  normals[3] = -1,0;
189 
190  for( size_t i = 0; i < normals.size(); ++i)
191  {
192  normals[i] /= sqrt( blitz::dot( normals[i], normals[i]));
193  }
194 
195  // make a rotation object from this contour
196  drawRotationObjectFromContour( origin, nx, ny, nz,
197  contour, normals, _nSegments);
198  }
199 
200 
201 
202  void drawArrow(const blitz::TinyVector<float,3>& origin,
203  const blitz::TinyVector<float,3>& dxyz)
204  {
205  // compute x- and y- unit vectors
206  float length = std::sqrt( blitz::dot( dxyz, dxyz));
207 
208  blitz::TinyVector<float,3> nx = dxyz / length;
209 
210  blitz::TinyVector<float,3> unitCol(0,0,1);
211  blitz::TinyVector<float,3> ny;
212  blitz::TinyVector<float,3> nz;
213 
214  ny = blitz::cross( nx, unitCol);
215  nz = blitz::cross( nx, ny);
216 
217  ny /= std::sqrt( blitz::dot( ny, ny));
218  nz /= std::sqrt( blitz::dot( nz, nz));
219 
220  // create the contour of the arrow
221  std::vector< blitz::TinyVector<float,2> > contour(5);
222  contour[0] = length, 0;
223  contour[1] = length-2*_headRadius, _headRadius;
224  contour[2] = length-2*_headRadius, _tailRadius;
225  contour[3] = 0, _tailRadius;
226  contour[4] = 0,0;
227 
228  std::vector< blitz::TinyVector<float,2> > normals(5);
229  normals[0] = 1, 0;
230  normals[1] = 1, 2;
231  normals[2] = -1, 0;
232  normals[3] = 0, 1;
233  normals[4] = -1,0;
234 
235  for( size_t i = 0; i < normals.size(); ++i)
236  {
237  normals[i] /= std::sqrt( blitz::dot( normals[i], normals[i]));
238  }
239 
240  // make a rotation object from this contour
241  drawRotationObjectFromContour( origin, nx, ny, nz,
242  contour, normals, _nSegments);
243  }
244 
245 
246  void drawEllipse(const blitz::TinyVector<float,3>& origin,
247  const blitz::TinyVector<float,3>& axis1,
248  const blitz::TinyVector<float,3>& axis2)
249  {
250  float radius1 = sqrt( blitz::dot( axis1, axis1));
251  float radius2 = sqrt( blitz::dot( axis2, axis2));
252 
253  blitz::TinyVector<float,3> ny = axis1 / radius1;
254  blitz::TinyVector<float,3> nz = axis2 / radius2;
255  blitz::TinyVector<float,3> nx;
256  nx = blitz::cross( ny, nz);
257 
258  // create the roational contour for the border of the disc
259  std::vector<blitz::TinyVector<float,2> > contour(_nSegments);
260  std::vector<blitz::TinyVector<float,2> > normals(_nSegments);
261  for( int i = 0; i < _nSegments; ++i)
262  {
263  float phi = float(i)/(_nSegments-1)*2*M_PI;
264  contour[i] =
265  _borderRadius*sin(phi),
266  _borderRadius*cos(phi);
267  //float phiNormal = (float(i)-0.5)/(_nSegments-1)*2*M_PI;
268  normals[i] = sin(phi), cos(phi);
269  }
270  // normalize normals
271  for( size_t i = 0; i < normals.size(); ++i)
272  {
273  normals[i] /= sqrt( blitz::dot( normals[i], normals[i]));
274  }
275 
276  // make a rotation object from this contour
277  for( int i = 0; i < _nEllipseSegments; ++i)
278  {
279  float phi1 = float(i)/_nEllipseSegments*2*M_PI;
280  float phi2 = float(i+1)/_nEllipseSegments*2*M_PI;
281  blitz::TinyVector<float,3> center1 =
282  radius1 * sin( phi1) * ny + radius2 * cos( phi1) * nz;
283 
284  blitz::TinyVector<float,3> rotNy1 =
285  center1 / sqrt( blitz::dot( center1, center1));
286 
287  blitz::TinyVector<float,3> center2 =
288  radius1 * sin( phi2) * ny + radius2 * cos( phi2) * nz;
289 
290  blitz::TinyVector<float,3> rotNy2 =
291  center2 / sqrt( blitz::dot( center2, center2));
292 
293  blitz::TinyVector<float,3> rotNyMean =
294  (rotNy1 + rotNy2)/2;
295  rotNyMean /= sqrt( blitz::dot( rotNyMean, rotNyMean));
296 
297  for( size_t j = 0; j < contour.size()-1; ++j)
298  {
299  blitz::TinyVector<float,3> p1 =
300  origin
301  + nx * contour[j](0)
302  + center1 + rotNy1 * contour[j](1);
303 
304  blitz::TinyVector<float,3> p2 =
305  origin
306  + nx * contour[j](0)
307  + center2 + rotNy2 * contour[j](1);
308 
309  blitz::TinyVector<float,3> p3 =
310  origin
311  + nx * contour[j+1](0)
312  + center1 + rotNy1 * contour[j+1](1);
313 
314  blitz::TinyVector<float,3> p4 =
315  origin
316  + nx * contour[j+1](0)
317  + center2 + rotNy2 * contour[j+1](1);
318 
319  blitz::TinyVector<float,3> normal =
320  nx * normals[j](0)
321  + rotNyMean * normals[j](1);
322 
323  // triangle 1 (p1 -- p2 -- p3)
324  if(blitz::any(p1 != p2))
325  {
326  drawTriangle( p1, p2, p3, normal);
327  }
328 
329  // triangle 2 (p2 -- p4 -- p3)
330  if(blitz::any(p3 != p4))
331  {
332  drawTriangle( p2, p4, p3, normal);
333  }
334 
335  }
336  }
337 
338 
339  }
340 
341  void drawSphere(const blitz::TinyVector<float,3>& center,
342  const float radius, const int halfCircleSampling = 10)
343  {
344  //sample points for one half circle
345  const int sampling = halfCircleSampling;
346  blitz::TinyVector<float,3> znormal(0,0,1);
347 
348  blitz::Array< blitz::TinyVector<float,3>, 1>
349  startHalfCircleSamples(sampling);
350  startHalfCircleSamples = blitz::TinyVector<float,3>(0,0,0);
351 
352  for (int i = 0; i < sampling; i++)
353  {
354  float factor = ( float(i)/float(sampling-1)) - 0.5;
355  float x = radius * cos( factor * M_PI);
356  float y = 0;
357  float z = radius * sin( factor * M_PI);
358  startHalfCircleSamples(i) = blitz::TinyVector<float,3>(x,y,z);
359  }
360 
361  blitz::Array< blitz::TinyVector<float,3>, 1>
362  halfCircleSamples(sampling);
363  halfCircleSamples = startHalfCircleSamples;
364  blitz::Array< blitz::TinyVector<float,3>, 1>
365  nextHalfCircleSamples(sampling);
366  nextHalfCircleSamples = blitz::TinyVector<float,3>(0,0,0);
367 
368  for (int i = 1; i < 2*sampling; ++i)
369  {
370  float factor = ( float(i)/float(2*sampling-1));
371  float phi = factor * 2*M_PI;
372  float cosPhi = cos(phi);
373  float sinPhi = sin(phi);
374 
375  // counter-clockwise rotation around z-axes
376  blitz::Array<float,2> rotMatrix(3,3);
377  rotMatrix(0,0) = cosPhi;
378  rotMatrix(0,1) = -sinPhi;
379  rotMatrix(0,2) = 0;
380 
381  rotMatrix(1,0) = sinPhi;
382  rotMatrix(1,1) = cosPhi;
383  rotMatrix(1,2) = 0;
384 
385  rotMatrix(2,0) = 0;
386  rotMatrix(2,1) = 0;
387  rotMatrix(2,2) = 1;
388 
389  for( int j = 0; j < startHalfCircleSamples.extent(0); ++j) {
390  nextHalfCircleSamples(j) =
391  mvMult(rotMatrix, startHalfCircleSamples(j));
392  }
393 
394  for ( int k = 0; k < sampling-2; ++k)
395  {
396  blitz::TinyVector<float,3> p1 = nextHalfCircleSamples(k);
397  blitz::TinyVector<float,3> p2 = nextHalfCircleSamples(k+1);
398  blitz::TinyVector<float,3> p3 = halfCircleSamples(k+1);
399  blitz::TinyVector<float,3> v1 = p2 - p1;
400  blitz::TinyVector<float,3> v2 = p3 - p1;
401  blitz::TinyVector<float,3> n = blitz::cross(v1,v2);
402  n /= sqrt(n(0)*n(0) + n(1)*n(1) + n(2)*n(2));
403  p1 += center; p2 += center; p3 += center;
404  drawTriangle(p1, p2, p3, n);
405  }
406 
407  for ( int l = 1; l < sampling-1; ++l)
408  {
409  blitz::TinyVector<float,3> p1 = halfCircleSamples(l);
410  blitz::TinyVector<float,3> p2 = nextHalfCircleSamples(l);
411  blitz::TinyVector<float,3> p3 = halfCircleSamples(l+1);
412  blitz::TinyVector<float,3> v1 = p2 - p1;
413  blitz::TinyVector<float,3> v2 = p3 - p1;
414  blitz::TinyVector<float,3> n = blitz::cross(v1,v2);
415  n /= sqrt(n(0)*n(0) + n(1)*n(1) + n(2)*n(2));
416  p1 += center; p2 += center; p3 += center;
417  drawTriangle(p1, p2, p3, n);
418  }
419 
420  halfCircleSamples = nextHalfCircleSamples;
421 
422  }
423 
424  }
425 
427  blitz::TinyVector<double,3> const &center,
428  blitz::Array<double,2> const &f)
429  {
430  for (ptrdiff_t thetaIdx = 0; thetaIdx < f.extent(0) - 1; ++thetaIdx)
431  {
432  double theta =
433  static_cast<double>(thetaIdx) /
434  static_cast<double>(f.extent(0) - 1) * M_PI;
435  double thetap =
436  static_cast<double>(thetaIdx + 1) /
437  static_cast<double>(f.extent(0) - 1) * M_PI;
438  for (ptrdiff_t phiIdx = 0; phiIdx < f.extent(1); ++phiIdx)
439  {
440  double phi = static_cast<double>(phiIdx) /
441  static_cast<double>(f.extent(1)) * 2.0 * M_PI;
442  double phip = (phiIdx == f.extent(1) - 1) ?
443  0.0 : (static_cast<double>(phiIdx + 1) /
444  static_cast<double>(f.extent(1)) * 2.0 * M_PI);
445 
446  double r = f(thetaIdx, phiIdx);
447  blitz::TinyVector<double,3> v1(
448  r * std::cos(theta),
449  r * std::sin(theta) * std::sin(phi),
450  r * std::sin(theta) * std::cos(phi));
451  v1 += center;
452 
453  r = f(thetaIdx + 1,
454  (phiIdx == f.extent(1) - 1) ? 0 : (phiIdx + 1));
455  blitz::TinyVector<double,3> v2(
456  r * std::cos(thetap),
457  r * std::sin(thetap) * std::sin(phip),
458  r * std::sin(thetap) * std::cos(phip));
459  v2 += center;
460 
461  r = f(thetaIdx + 1, phiIdx);
462  blitz::TinyVector<double,3> v3(
463  r * std::cos(thetap),
464  r * std::sin(thetap) * std::sin(phi),
465  r * std::sin(thetap) * std::cos(phi));
466  v3 += center;
467 
468  r = f(thetaIdx,
469  (phiIdx == f.extent(1) - 1) ? 0 : (phiIdx + 1));
470  blitz::TinyVector<double,3> v4(
471  r * std::cos(theta),
472  r * std::sin(theta) * std::sin(phip),
473  r * std::sin(theta) * std::cos(phip));
474  v4 += center;
475 
476  blitz::TinyVector<double,3> d21(v2 - v1);
477  blitz::TinyVector<double,3> d31(v3 - v1);
478  blitz::TinyVector<double,3> d41(v4 - v1);
479 
480  blitz::TinyVector<double,3> n(
481  blitz::cross(
482  d21 / std::sqrt(blitz::dot(d21, d21)),
483  d31 / std::sqrt(blitz::dot(d31, d31))));
484 
485  drawTriangle(v1, v2, v3, n);
486 
487  n = blitz::cross(
488  d41 / std::sqrt(blitz::dot(d41, d41)),
489  d21 / std::sqrt(blitz::dot(d21, d21))));
490 
491  drawTriangle(v1, v4, v2, n);
492  }
493  }
494  }
495 
496  };
497 
498 }
499 
500 #endif
void drawArrow(const blitz::TinyVector< float, 3 > &origin, const blitz::TinyVector< float, 3 > &dxyz)
void drawTriangle(const blitz::TinyVector< float, 3 > &p1, const blitz::TinyVector< float, 3 > &p2, const blitz::TinyVector< float, 3 > &p3, const blitz::TinyVector< float, 3 > &normal=blitz::TinyVector< float, 3 >(0, 0, 0))
void drawEllipse(const blitz::TinyVector< float, 3 > &origin, const blitz::TinyVector< float, 3 > &axis1, const blitz::TinyVector< float, 3 > &axis2)
bool any(blitz::TinyMatrix< bool, NRows, NColumns > const &matrix)
any() reduction for boolean blitz::TinyMatrix.
void drawSphere(const blitz::TinyVector< float, 3 > &center, const float radius, const int halfCircleSampling=10)
void drawRotationObjectFromContour(const blitz::TinyVector< float, 3 > &origin, const blitz::TinyVector< float, 3 > &nx, const blitz::TinyVector< float, 3 > &ny, const blitz::TinyVector< float, 3 > &nz, const std::vector< blitz::TinyVector< float, 2 > > &contour, const std::vector< blitz::TinyVector< float, 2 > > &normals, const int nAngles)
blitz::TinyVector< BaseT, Dim > mvMult(MatrixT const &m, blitz::TinyVector< BaseT, Dim > const &v)
Quadratic Matrix * Vector product (generic).
void drawStarshapedSurface(blitz::TinyVector< double, 3 > const &center, blitz::Array< double, 2 > const &f)
STLFileWriter(const std::string &fileName, const std::string &objectName="dummy")
void drawCylinder(const blitz::TinyVector< float, 3 > &origin, const blitz::TinyVector< float, 3 > &dxyz, double radius=-1)
double dot(blitz::Array< DataT, 1 > const &vec1, blitz::Array< DataT, 1 > const &vec2)
Calculation of the inner product (dot product) of the given vectors.