gmsh-TingyuanDoc  0.1
An Open-Source Timing-driven Analytical Mixed-size FPGA Placer
SOrientedBoundingBox.cpp
Go to the documentation of this file.
1 // Gmsh - Copyright (C) 1997-2022 C. Geuzaine, J.-F. Remacle
2 //
3 // See the LICENSE.txt file in the Gmsh root directory for license information.
4 // Please report all issues on https://gitlab.onelab.info/gmsh/gmsh/issues.
5 //
6 // Contributor(s):
7 // Bastien Gorissen
8 //
9 
10 #include <algorithm>
11 #include <cmath>
12 #include <time.h>
13 #include <limits>
14 #include "GmshConfig.h"
15 #include "SOrientedBoundingBox.h"
16 #include "fullMatrix.h"
17 #include "SBoundingBox3d.h"
18 #include "fullMatrix.h"
19 
20 #if defined(HAVE_MESH)
21 #include "DivideAndConquer.h"
22 #endif
23 
25  : center({{0.0, 0.0}}), size({{0.0, 0.0}}), axisX({{0.0, 0.0}}),
26  axisY({{0.0, 0.0}})
27 {
28 }
29 
30 double SOrientedBoundingRectangle::area() { return size[0] * size[1]; }
31 
33 {
34  double dx = 0.5 * size[0];
35  double dy = 0.5 * size[1];
36  double dz = 0.5 * size[2];
37 
38  p1x = center[0] - (axisX[0] * dx) - (axisY[0] * dy) - (axisZ[0] * dz);
39  p1y = center[1] - (axisX[1] * dx) - (axisY[1] * dy) - (axisZ[1] * dz);
40  p1z = center[2] - (axisX[2] * dx) - (axisY[2] * dy) - (axisZ[2] * dz);
41 
42  p2x = center[0] + (axisX[0] * dx) - (axisY[0] * dy) - (axisZ[0] * dz);
43  p2y = center[1] + (axisX[1] * dx) - (axisY[1] * dy) - (axisZ[1] * dz);
44  p2z = center[2] + (axisX[2] * dx) - (axisY[2] * dy) - (axisZ[2] * dz);
45 
46  p3x = center[0] - (axisX[0] * dx) + (axisY[0] * dy) - (axisZ[0] * dz);
47  p3y = center[1] - (axisX[1] * dx) + (axisY[1] * dy) - (axisZ[1] * dz);
48  p3z = center[2] - (axisX[2] * dx) + (axisY[2] * dy) - (axisZ[2] * dz);
49 
50  p4x = center[0] + (axisX[0] * dx) + (axisY[0] * dy) - (axisZ[0] * dz);
51  p4y = center[1] + (axisX[1] * dx) + (axisY[1] * dy) - (axisZ[1] * dz);
52  p4z = center[2] + (axisX[2] * dx) + (axisY[2] * dy) - (axisZ[2] * dz);
53 
54  p5x = center[0] - (axisX[0] * dx) - (axisY[0] * dy) + (axisZ[0] * dz);
55  p5y = center[1] - (axisX[1] * dx) - (axisY[1] * dy) + (axisZ[1] * dz);
56  p5z = center[2] - (axisX[2] * dx) - (axisY[2] * dy) + (axisZ[2] * dz);
57 
58  p6x = center[0] + (axisX[0] * dx) - (axisY[0] * dy) + (axisZ[0] * dz);
59  p6y = center[1] + (axisX[1] * dx) - (axisY[1] * dy) + (axisZ[1] * dz);
60  p6z = center[2] + (axisX[2] * dx) - (axisY[2] * dy) + (axisZ[2] * dz);
61 
62  p7x = center[0] - (axisX[0] * dx) + (axisY[0] * dy) + (axisZ[0] * dz);
63  p7y = center[1] - (axisX[1] * dx) + (axisY[1] * dy) + (axisZ[1] * dz);
64  p7z = center[2] - (axisX[2] * dx) + (axisY[2] * dy) + (axisZ[2] * dz);
65 
66  p8x = center[0] + (axisX[0] * dx) + (axisY[0] * dy) + (axisZ[0] * dz);
67  p8y = center[1] + (axisX[1] * dx) + (axisY[1] * dy) + (axisZ[1] * dz);
68  p8z = center[2] + (axisX[2] * dx) + (axisY[2] * dy) + (axisZ[2] * dz);
69 }
70 
72 {
73  center = SVector3();
74  size = SVector3();
75  axisX = SVector3();
76  axisY = SVector3();
77  axisZ = SVector3();
78  fillp();
79 }
80 
82  double sizeY, double sizeZ,
83  const SVector3 &axisX_,
84  const SVector3 &axisY_,
85  const SVector3 &axisZ_)
86 {
87  center = center_;
88  size = SVector3(sizeX, sizeY, sizeZ);
89  axisX = axisX_;
90  axisX.normalize();
91  axisY = axisY_;
92  axisY.normalize();
93  axisZ = axisZ_;
94  axisZ.normalize();
95  fillp();
96 }
97 
99 {
100  size = other->getSize();
101  axisX = other->getAxis(0);
102  axisY = other->getAxis(1);
103  axisZ = other->getAxis(2);
104  center = other->getCenter();
105  fillp();
106 }
107 
109 {
110  SVector3 ret;
111  switch(axis) {
112  case 0: ret = axisX; break;
113  case 1: ret = axisY; break;
114  case 2: ret = axisZ; break;
115  }
116  return ret;
117 }
118 
120 {
121  SVector3 collide_axes[15];
122  for(int i = 0; i < 3; i++) {
123  collide_axes[i] = getAxis(i);
124  collide_axes[i + 3] = obb.getAxis(i);
125  }
126 
127  SVector3 sizes[2];
128  sizes[0] = getSize();
129  sizes[1] = obb.getSize();
130 
131  for(std::size_t i = 0; i < 3; i++) {
132  for(std::size_t j = 3; j < 6; j++) {
133  collide_axes[3 * i + j + 3] = crossprod(collide_axes[i], collide_axes[j]);
134  }
135  }
136  SVector3 T = obb.getCenter() - getCenter();
137 
138  for(std::size_t i = 0; i < 15; i++) {
139  double val = 0.0;
140  for(std::size_t j = 0; j < 6; j++) {
141  val += 0.5 * (sizes[j < 3 ? 0 : 1])(j % 3) *
142  std::abs(dot(collide_axes[j], collide_axes[i]));
143  }
144  if(std::abs(dot(collide_axes[i], T)) > val) { return false; }
145  }
146  return true;
147 }
148 
150 SOrientedBoundingBox::buildOBB(std::vector<SPoint3> &vertices)
151 {
152 #if defined(HAVE_MESH)
153 
154  int num_vertices = vertices.size();
155  // First organize the data
156 
157  std::set<SPoint3> unique;
158  unique.insert(vertices.begin(), vertices.end());
159 
160  num_vertices = unique.size();
161  fullMatrix<double> data(3, num_vertices);
162 
163  fullVector<double> mean(3);
164  fullVector<double> vmins(3);
165  fullVector<double> vmaxs(3);
166 
167  mean.setAll(0);
168  vmins.setAll(DBL_MAX);
169  vmaxs.setAll(-DBL_MAX);
170 
171  size_t idx = 0;
172  for(auto uIter = unique.begin(); uIter != unique.end(); ++uIter) {
173  const SPoint3 &pp = *uIter;
174  for(int d = 0; d < 3; d++) {
175  data(d, idx) = pp[d];
176  vmins(d) = std::min(vmins(d), pp[d]);
177  vmaxs(d) = std::max(vmaxs(d), pp[d]);
178  mean(d) += pp[d];
179  }
180  idx++;
181  }
182 
183  for(int i = 0; i < 3; i++) { mean(i) /= num_vertices; }
184 
185  // Get the deviation from the mean
186  fullMatrix<double> B(3, num_vertices);
187  for(int i = 0; i < 3; i++) {
188  for(int j = 0; j < num_vertices; j++) { B(i, j) = data(i, j) - mean(i); }
189  }
190 
191  // Compute the covariance matrix
192  fullMatrix<double> covariance(3, 3);
193  B.mult(B.transpose(), covariance);
194  covariance.scale(1. / (num_vertices - 1));
195  /*
196  Msg::Debug("Covariance matrix");
197  Msg::Debug("%f %f %f", covariance(0,0),covariance(0,1),covariance(0,2) );
198  Msg::Debug("%f %f %f", covariance(1,0),covariance(1,1),covariance(1,2) );
199  Msg::Debug("%f %f %f", covariance(2,0),covariance(2,1),covariance(2,2) );
200  */
201  for(int i = 0; i < 3; i++) {
202  for(int j = 0; j < 3; j++) {
203  if(std::abs(covariance(i, j)) < 10e-16) covariance(i, j) = 0;
204  }
205  }
206 
207  fullMatrix<double> left_eigv(3, 3);
208  fullMatrix<double> right_eigv(3, 3);
209  fullVector<double> real_eig(3);
210  fullVector<double> img_eig(3);
211  covariance.eig(real_eig, img_eig, left_eigv, right_eigv, true);
212 
213  // Now, project the data in the new basis.
214  fullMatrix<double> projected(3, num_vertices);
215  left_eigv.transpose().mult(data, projected);
216  // Get the size of the box in the new direction
217  fullVector<double> mins(3);
218  fullVector<double> maxs(3);
219  for(int i = 0; i < 3; i++) {
220  mins(i) = DBL_MAX;
221  maxs(i) = -DBL_MAX;
222  for(int j = 0; j < num_vertices; j++) {
223  maxs(i) = std::max(maxs(i), projected(i, j));
224  mins(i) = std::min(mins(i), projected(i, j));
225  }
226  }
227 
228  // double means[3];
229  double sizes[3];
230 
231  // Note: the size is computed in the box's coordinates!
232  for(int i = 0; i < 3; i++) {
233  sizes[i] = maxs(i) - mins(i);
234  // means[i] = (maxs(i) - mins(i)) / 2.;
235  }
236 
237  if(sizes[0] == 0 && sizes[1] == 0) {
238  // Entity is a straight line...
240  SVector3 Axis1;
241  SVector3 Axis2;
242  SVector3 Axis3;
243 
244  Axis1[0] = left_eigv(0, 0);
245  Axis1[1] = left_eigv(1, 0);
246  Axis1[2] = left_eigv(2, 0);
247  Axis2[0] = left_eigv(0, 1);
248  Axis2[1] = left_eigv(1, 1);
249  Axis2[2] = left_eigv(2, 1);
250  Axis3[0] = left_eigv(0, 2);
251  Axis3[1] = left_eigv(1, 2);
252  Axis3[2] = left_eigv(2, 2);
253 
254  center[0] = (vmaxs(0) + vmins(0)) / 2.0;
255  center[1] = (vmaxs(1) + vmins(1)) / 2.0;
256  center[2] = (vmaxs(2) + vmins(2)) / 2.0;
257 
258  return new SOrientedBoundingBox(center, sizes[0], sizes[1], sizes[2], Axis1,
259  Axis2, Axis3);
260  }
261 
262  // We take the smallest component, then project the data on the plane defined
263  // by the other twos
264 
265  int smallest_comp = 0;
266  if(sizes[0] <= sizes[1] && sizes[0] <= sizes[2])
267  smallest_comp = 0;
268  else if(sizes[1] <= sizes[0] && sizes[1] <= sizes[2])
269  smallest_comp = 1;
270  else if(sizes[2] <= sizes[0] && sizes[2] <= sizes[1])
271  smallest_comp = 2;
272 
273  // The projection has been done circa line 161.
274  // We just ignore the coordinate corresponding to smallest_comp.
275  std::vector<SPoint2 *> points;
276  for(int i = 0; i < num_vertices; i++) {
277  SPoint2 *p = new SPoint2(projected(smallest_comp == 0 ? 1 : 0, i),
278  projected(smallest_comp == 2 ? 1 : 2, i));
279  bool keep = true;
280  for(auto point = points.begin(); point != points.end(); point++) {
281  if(std::abs((*p)[0] - (**point)[0]) < 10e-10 &&
282  std::abs((*p)[1] - (**point)[1]) < 10e-10) {
283  keep = false;
284  break;
285  }
286  }
287  if(keep) { points.push_back(p); }
288  else {
289  delete p;
290  }
291  }
292 
293  // Find the convex hull from a delaunay triangulation of the points
294  DocRecord record(points.size());
295  record.numPoints = points.size();
296  srand((unsigned)time(nullptr));
297  for(std::size_t i = 0; i < points.size(); i++) {
298  record.points[i].where.h =
299  points[i]->x() + (10e-6) * sizes[smallest_comp == 0 ? 1 : 0] *
300  (-0.5 + ((double)rand()) / RAND_MAX);
301  record.points[i].where.v =
302  points[i]->y() + (10e-6) * sizes[smallest_comp == 2 ? 1 : 0] *
303  (-0.5 + ((double)rand()) / RAND_MAX);
304  record.points[i].adjacent = nullptr;
305  }
306 
307  try {
308  record.MakeMeshWithPoints();
309  } catch(std::runtime_error &e) {
310  Msg::Error("%s", e.what());
311  }
312 
313  std::vector<Segment> convex_hull;
314  for(int i = 0; i < record.numTriangles; i++) {
315  Segment segs[3];
316  segs[0].from = record.triangles[i].a;
317  segs[0].to = record.triangles[i].b;
318  segs[1].from = record.triangles[i].b;
319  segs[1].to = record.triangles[i].c;
320  segs[2].from = record.triangles[i].c;
321  segs[2].to = record.triangles[i].a;
322 
323  for(int j = 0; j < 3; j++) {
324  bool okay = true;
325  for(auto seg = convex_hull.begin(); seg != convex_hull.end(); seg++) {
326  if(((*seg).from == segs[j].from && (*seg).from == segs[j].to)
327  // FIXME:
328  // || ((*seg).from == segs[j].to && (*seg).from == segs[j].from)
329  ) {
330  convex_hull.erase(seg);
331  okay = false;
332  break;
333  }
334  }
335  if(okay) { convex_hull.push_back(segs[j]); }
336  }
337  }
338 
339  // Now, examinate all the directions given by the edges of the convex hull
340  // to find the one that lets us build the least-area bounding rectangle for
341  // then points.
342  fullVector<double> axis(2);
343  axis(0) = 1;
344  axis(1) = 0;
345  fullVector<double> axis2(2);
346  axis2(0) = 0;
347  axis2(1) = 1;
348  SOrientedBoundingRectangle least_rectangle;
349  least_rectangle.center[0] = 0.0;
350  least_rectangle.center[1] = 0.0;
351  least_rectangle.size[0] = -1.0;
352  least_rectangle.size[1] = 1.0;
353 
354  fullVector<double> segment(2);
355  fullMatrix<double> rotation(2, 2);
356 
357  for(auto seg = convex_hull.begin(); seg != convex_hull.end(); seg++) {
358  // segment(0) = record.points[(*seg).from].where.h -
359  // record.points[(*seg).to].where.h; segment(1) =
360  // record.points[(*seg).from].where.v - record.points[(*seg).to].where.v;
361  segment(0) = points[(*seg).from]->x() - points[(*seg).to]->x();
362  segment(1) = points[(*seg).from]->y() - points[(*seg).to]->y();
363  segment.scale(1.0 / segment.norm());
364 
365  double cosine = axis(0) * segment(0) + segment(1) * axis(1);
366  double sine = axis(1) * segment(0) - segment(1) * axis(0);
367  // double sine = axis(0)*segment(1) - segment(0)*axis(1);
368 
369  rotation(0, 0) = cosine;
370  rotation(0, 1) = sine;
371  rotation(1, 0) = -sine;
372  rotation(1, 1) = cosine;
373 
374  auto max_x = -std::numeric_limits<double>::max();
375  auto min_x = std::numeric_limits<double>::max();
376  auto max_y = -std::numeric_limits<double>::max();
377  auto min_y = std::numeric_limits<double>::max();
378 
379  for(int i = 0; i < record.numPoints; i++) {
380  fullVector<double> pnt(2);
381  // pnt(0) = record.points[i].where.h;
382  // pnt(1) = record.points[i].where.v;
383  pnt(0) = points[i]->x();
384  pnt(1) = points[i]->y();
385 
386  fullVector<double> rot_pnt(2);
387  rotation.mult(pnt, rot_pnt);
388 
389  if(rot_pnt(0) < min_x) min_x = rot_pnt(0);
390  if(rot_pnt(0) > max_x) max_x = rot_pnt(0);
391  if(rot_pnt(1) < min_y) min_y = rot_pnt(1);
392  if(rot_pnt(1) > max_y) max_y = rot_pnt(1);
393  }
394 
395 
396  fullVector<double> center_rot(2);
397  fullVector<double> center_before_rot(2);
398  center_before_rot(0) = (max_x + min_x) / 2.0;
399  center_before_rot(1) = (max_y + min_y) / 2.0;
400  fullMatrix<double> rotation_inv(2, 2);
401 
402  rotation_inv(0, 0) = cosine;
403  rotation_inv(0, 1) = -sine;
404  rotation_inv(1, 0) = sine;
405  rotation_inv(1, 1) = cosine;
406 
407  rotation_inv.mult(center_before_rot, center_rot);
408 
409  fullVector<double> axis_rot1(2);
410  fullVector<double> axis_rot2(2);
411 
412  rotation_inv.mult(axis, axis_rot1);
413  rotation_inv.mult(axis2, axis_rot2);
414 
415  if((least_rectangle.area() == -1) ||
416  (max_x - min_x) * (max_y - min_y) < least_rectangle.area()) {
417  least_rectangle.size[0] = max_x - min_x;
418  least_rectangle.size[1] = max_y - min_y;
419  least_rectangle.center[0] = (max_x + min_x) / 2.0;
420  least_rectangle.center[1] = (max_y + min_y) / 2.0;
421  least_rectangle.center[0] = center_rot(0);
422  least_rectangle.center[1] = center_rot(1);
423  least_rectangle.axisX[0] = axis_rot1(0);
424  least_rectangle.axisX[1] = axis_rot1(1);
425  // least_rectangle.axisX[0] = segment(0);
426  // least_rectangle.axisX[1] = segment(1);
427  least_rectangle.axisY[0] = axis_rot2(0);
428  least_rectangle.axisY[1] = axis_rot2(1);
429  }
430  }
431 
432  auto min_pca = std::numeric_limits<double>::max();
433  auto max_pca = -std::numeric_limits<double>::max();
434  for(int i = 0; i < num_vertices; i++) {
435  min_pca = std::min(min_pca, projected(smallest_comp, i));
436  max_pca = std::max(max_pca, projected(smallest_comp, i));
437  }
438  double center_pca = (max_pca + min_pca) / 2.0;
439  double size_pca = (max_pca - min_pca);
440 
441  double raw_data[3][5];
442  raw_data[0][0] = size_pca;
443  raw_data[1][0] = least_rectangle.size[0];
444  raw_data[2][0] = least_rectangle.size[1];
445 
446  raw_data[0][1] = center_pca;
447  raw_data[1][1] = least_rectangle.center[0];
448  raw_data[2][1] = least_rectangle.center[1];
449 
450  for(int i = 0; i < 3; i++) {
451  raw_data[0][2 + i] = left_eigv(i, smallest_comp);
452  raw_data[1][2 + i] =
453  least_rectangle.axisX[0] * left_eigv(i, smallest_comp == 0 ? 1 : 0) +
454  least_rectangle.axisX[1] * left_eigv(i, smallest_comp == 2 ? 1 : 2);
455  raw_data[2][2 + i] =
456  least_rectangle.axisY[0] * left_eigv(i, smallest_comp == 0 ? 1 : 0) +
457  least_rectangle.axisY[1] * left_eigv(i, smallest_comp == 2 ? 1 : 2);
458  }
459  // Msg::Info("Test 1 : %f
460  // %f",least_rectangle.center[0],least_rectangle.center[1]);
461  // Msg::Info("Test 2 : %f
462  // %f",least_rectangle.axisY[0],least_rectangle.axisY[1]);
463 
464  int tri[3];
465 
466  if(size_pca > least_rectangle.size[0]) {
467  // P > R0
468  if(size_pca > least_rectangle.size[1]) {
469  // P > R1
470  tri[0] = 0;
471  if(least_rectangle.size[0] > least_rectangle.size[1]) {
472  // R0 > R1
473  tri[1] = 1;
474  tri[2] = 2;
475  }
476  else {
477  // R1 > R0
478  tri[1] = 2;
479  tri[2] = 1;
480  }
481  }
482  else {
483  // P < R1
484  tri[0] = 2;
485  tri[1] = 0;
486  tri[2] = 1;
487  }
488  }
489  else { // P < R0
490  if(size_pca < least_rectangle.size[1]) {
491  // P < R1
492  tri[2] = 0;
493  if(least_rectangle.size[0] > least_rectangle.size[1]) {
494  tri[0] = 1;
495  tri[1] = 2;
496  }
497  else {
498  tri[0] = 2;
499  tri[1] = 1;
500  }
501  }
502  else {
503  tri[0] = 1;
504  tri[1] = 0;
505  tri[2] = 2;
506  }
507  }
508 
509  SVector3 size;
511  SVector3 Axis1;
512  SVector3 Axis2;
513  SVector3 Axis3;
514 
515  for(int i = 0; i < 3; i++) {
516  size[i] = raw_data[tri[i]][0];
517  center[i] = raw_data[tri[i]][1];
518  Axis1[i] = raw_data[tri[0]][2 + i];
519  Axis2[i] = raw_data[tri[1]][2 + i];
520  Axis3[i] = raw_data[tri[2]][2 + i];
521  }
522 
523  SVector3 aux1;
524  SVector3 aux2;
525  SVector3 aux3;
526  for(int i = 0; i < 3; i++) {
527  aux1(i) = left_eigv(i, smallest_comp);
528  aux2(i) = left_eigv(i, smallest_comp == 0 ? 1 : 0);
529  aux3(i) = left_eigv(i, smallest_comp == 2 ? 1 : 2);
530  }
531  center = aux1 * center_pca + aux2 * least_rectangle.center[0] +
532  aux3 * least_rectangle.center[1];
533  // center[1] = -center[1];
534 
535  /*
536  Msg::Info("Box center : %f %f %f",center[0],center[1],center[2]);
537  Msg::Info("Box size : %f %f %f",size[0],size[1],size[2]);
538  Msg::Info("Box axis 1 : %f %f %f",Axis1[0],Axis1[1],Axis1[2]);
539  Msg::Info("Box axis 2 : %f %f %f",Axis2[0],Axis2[1],Axis2[2]);
540  Msg::Info("Box axis 3 : %f %f %f",Axis3[0],Axis3[1],Axis3[2]);
541 
542  Msg::Info("Volume : %f", size[0]*size[1]*size[2]);
543  */
544 
545  return new SOrientedBoundingBox(center, size[0], size[1], size[2], Axis1,
546  Axis2, Axis3);
547 #else
548  Msg::Error("SOrientedBoundingBox requires mesh module");
549  return 0;
550 #endif
551 }
552 
554  SOrientedBoundingBox &obb2)
555 {
556  // "center term"
557  double center_term = norm(obb1.getCenter() - obb2.getCenter());
558 
559  // "size term"
560  double size_term = 0.0;
561  for(int i = 0; i < 3; i++) {
562  if(obb1.getSize()(i) + obb2.getSize()(i) != 0) {
563  size_term += std::abs(obb1.getSize()(i) - obb2.getSize()(i)) /
564  (obb1.getSize()(i) + obb2.getSize()(i));
565  }
566  }
567 
568  // "orientation term"
569  double orientation_term = 0.0;
570  for(int i = 0; i < 3; i++) {
571  orientation_term += 1 - std::abs(dot(obb1.getAxis(i), obb2.getAxis(i)));
572  }
573 
574  return center_term + size_term + orientation_term;
575 }
crossprod
SVector3 crossprod(const SVector3 &a, const SVector3 &b)
Definition: SVector3.h:150
SOrientedBoundingBox::axisY
SVector3 axisY
Definition: SOrientedBoundingBox.h:89
SOrientedBoundingBox::p3y
double p3y
Definition: SOrientedBoundingBox.h:75
dot
double dot(const SVector3 &a, const SMetric3 &m, const SVector3 &b)
Definition: STensor3.h:71
DocRecord::numPoints
int numPoints
Definition: DivideAndConquer.h:93
SOrientedBoundingBox::size
SVector3 size
Definition: SOrientedBoundingBox.h:87
SOrientedBoundingBox::p7y
double p7y
Definition: SOrientedBoundingBox.h:79
fullVector< double >
SPoint2
Definition: SPoint2.h:12
SOrientedBoundingBox::getSize
SVector3 getSize() const
Definition: SOrientedBoundingBox.h:51
Msg::Error
static void Error(const char *fmt,...)
Definition: GmshMessage.cpp:482
fullVector::norm
scalar norm() const
SOrientedBoundingBox::p5x
double p5x
Definition: SOrientedBoundingBox.h:77
SOrientedBoundingRectangle::center
std::array< double, 2 > center
Definition: SOrientedBoundingBox.h:27
PointRecord::where
DPoint where
Definition: DivideAndConquer.h:28
SPoint3
Definition: SPoint3.h:14
SOrientedBoundingBox::axisZ
SVector3 axisZ
Definition: SOrientedBoundingBox.h:90
fullMatrix::scale
void scale(const scalar s)
Definition: fullMatrix.h:447
SOrientedBoundingBox.h
SVector3
Definition: SVector3.h:16
SOrientedBoundingBox
Definition: SOrientedBoundingBox.h:33
fullVector::scale
void scale(const scalar s)
Definition: fullMatrix.h:144
fullMatrix::transpose
fullMatrix< scalar > transpose() const
Definition: fullMatrix.h:558
DocRecord::points
PointRecord * points
Definition: DivideAndConquer.h:95
Segment::to
PointNumero to
Definition: DivideAndConquer.h:57
DivideAndConquer.h
SOrientedBoundingRectangle::size
std::array< double, 2 > size
Definition: SOrientedBoundingBox.h:28
DocRecord
Definition: DivideAndConquer.h:64
SOrientedBoundingBox::p6z
double p6z
Definition: SOrientedBoundingBox.h:78
SOrientedBoundingBox::p8x
double p8x
Definition: SOrientedBoundingBox.h:80
SOrientedBoundingBox::p8y
double p8y
Definition: SOrientedBoundingBox.h:80
SOrientedBoundingBox::p4z
double p4z
Definition: SOrientedBoundingBox.h:76
SOrientedBoundingBox::p4x
double p4x
Definition: SOrientedBoundingBox.h:76
SOrientedBoundingBox::fillp
void fillp()
Definition: SOrientedBoundingBox.cpp:32
Triangle::c
PointNumero c
Definition: DivideAndConquer.h:61
PointRecord::adjacent
DListPeek adjacent
Definition: DivideAndConquer.h:29
SOrientedBoundingBox::p5z
double p5z
Definition: SOrientedBoundingBox.h:77
fullMatrix< double >
DocRecord::numTriangles
int numTriangles
Definition: DivideAndConquer.h:96
Triangle::a
PointNumero a
Definition: DivideAndConquer.h:61
SOrientedBoundingRectangle::area
double area()
Definition: SOrientedBoundingBox.cpp:30
SOrientedBoundingBox::p5y
double p5y
Definition: SOrientedBoundingBox.h:77
SOrientedBoundingRectangle::axisX
std::array< double, 2 > axisX
Definition: SOrientedBoundingBox.h:29
SOrientedBoundingBox::center
SVector3 center
Definition: SOrientedBoundingBox.h:86
SBoundingBox3d.h
norm
void norm(const double *vec, double *norm)
Definition: gmshLevelset.cpp:202
SOrientedBoundingBox::p2y
double p2y
Definition: SOrientedBoundingBox.h:74
SOrientedBoundingBox::intersects
bool intersects(SOrientedBoundingBox &obb) const
Definition: SOrientedBoundingBox.cpp:119
DPoint::h
double h
Definition: DivideAndConquer.h:24
SOrientedBoundingBox::getCenter
const SVector3 & getCenter() const
Definition: SOrientedBoundingBox.h:47
SOrientedBoundingRectangle
Definition: SOrientedBoundingBox.h:20
SOrientedBoundingBox::p1z
double p1z
Definition: SOrientedBoundingBox.h:73
fullMatrix::eig
bool eig(fullVector< double > &eigenValReal, fullVector< double > &eigenValImag, fullMatrix< scalar > &leftEigenVect, fullMatrix< scalar > &rightEigenVect, bool sortRealPart=false)
Definition: fullMatrix.h:727
SOrientedBoundingBox::SOrientedBoundingBox
SOrientedBoundingBox()
Definition: SOrientedBoundingBox.cpp:71
SOrientedBoundingBox::p4y
double p4y
Definition: SOrientedBoundingBox.h:76
Segment
Definition: DivideAndConquer.h:55
SOrientedBoundingRectangle::SOrientedBoundingRectangle
SOrientedBoundingRectangle()
Definition: SOrientedBoundingBox.cpp:24
Segment::from
PointNumero from
Definition: DivideAndConquer.h:56
SOrientedBoundingRectangle::axisY
std::array< double, 2 > axisY
Definition: SOrientedBoundingBox.h:30
DocRecord::triangles
Triangle * triangles
Definition: DivideAndConquer.h:97
DPoint::v
double v
Definition: DivideAndConquer.h:23
SOrientedBoundingBox::p6y
double p6y
Definition: SOrientedBoundingBox.h:78
SOrientedBoundingBox::compare
static double compare(SOrientedBoundingBox &obb1, SOrientedBoundingBox &obb2)
Definition: SOrientedBoundingBox.cpp:553
SOrientedBoundingBox::p8z
double p8z
Definition: SOrientedBoundingBox.h:80
SOrientedBoundingBox::axisX
SVector3 axisX
Definition: SOrientedBoundingBox.h:88
SOrientedBoundingBox::buildOBB
static SOrientedBoundingBox * buildOBB(std::vector< SPoint3 > &vertices)
Definition: SOrientedBoundingBox.cpp:150
SOrientedBoundingBox::p7z
double p7z
Definition: SOrientedBoundingBox.h:79
fullMatrix::mult
void mult(const fullVector< scalar > &x, fullVector< scalar > &y) const
Definition: fullMatrix.h:487
SOrientedBoundingBox::p3z
double p3z
Definition: SOrientedBoundingBox.h:75
SOrientedBoundingBox::p2x
double p2x
Definition: SOrientedBoundingBox.h:74
Triangle::b
PointNumero b
Definition: DivideAndConquer.h:61
SOrientedBoundingBox::p7x
double p7x
Definition: SOrientedBoundingBox.h:79
SOrientedBoundingBox::p6x
double p6x
Definition: SOrientedBoundingBox.h:78
SOrientedBoundingBox::p1x
double p1x
Definition: SOrientedBoundingBox.h:73
SOrientedBoundingBox::p3x
double p3x
Definition: SOrientedBoundingBox.h:75
DocRecord::MakeMeshWithPoints
void MakeMeshWithPoints()
Definition: DivideAndConquer.cpp:855
fullVector::setAll
void setAll(const scalar &m)
Definition: fullMatrix.h:158
fullMatrix.h
SOrientedBoundingBox::p1y
double p1y
Definition: SOrientedBoundingBox.h:73
SOrientedBoundingBox::getAxis
SVector3 getAxis(int axis) const
valid values for axis are 0 (X-axis), 1 (Y-axis) or 2 (Z-axis)
Definition: SOrientedBoundingBox.cpp:108
SOrientedBoundingBox::p2z
double p2z
Definition: SOrientedBoundingBox.h:74
SVector3::normalize
double normalize()
Definition: SVector3.h:38