VTK  9.1.0
vtkBoundingBox.h
Go to the documentation of this file.
1/*=========================================================================
2
3Program: Visualization Toolkit
4Module: vtkBoundingBox.h
5
6Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7All rights reserved.
8See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10This software is distributed WITHOUT ANY WARRANTY; without even
11the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
37#ifndef vtkBoundingBox_h
38#define vtkBoundingBox_h
39#include "vtkCommonDataModelModule.h" // For export macro
40#include "vtkSystemIncludes.h"
41#include <atomic> // For threaded bounding box computation
42
43class vtkPoints;
44
45class VTKCOMMONDATAMODEL_EXPORT vtkBoundingBox
46{
47public:
49
54 vtkBoundingBox(const double bounds[6]);
55 vtkBoundingBox(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
57
61 vtkBoundingBox(const vtkBoundingBox& bbox);
62
66 vtkBoundingBox& operator=(const vtkBoundingBox& bbox);
67
69
72 bool operator==(const vtkBoundingBox& bbox) const;
73 bool operator!=(const vtkBoundingBox& bbox) const;
75
77
81 void SetBounds(const double bounds[6]);
82 void SetBounds(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
84
86
93 static void ComputeBounds(vtkPoints* pts, double bounds[6]);
94 static void ComputeBounds(vtkPoints* pts, const unsigned char* ptUses, double bounds[6]);
95 static void ComputeBounds(
96 vtkPoints* pts, const std::atomic<unsigned char>* ptUses, double bounds[6]);
97 void ComputeBounds(vtkPoints* pts) { this->ComputeBounds(pts, (unsigned char*)nullptr); }
98 void ComputeBounds(vtkPoints* pts, unsigned char* ptUses)
99 {
100 double bds[6];
101 vtkBoundingBox::ComputeBounds(pts, ptUses, bds);
102 this->MinPnt[0] = bds[0];
103 this->MinPnt[1] = bds[2];
104 this->MinPnt[2] = bds[4];
105 this->MaxPnt[0] = bds[1];
106 this->MaxPnt[1] = bds[3];
107 this->MaxPnt[2] = bds[5];
108 }
110
112
117 vtkPoints* points, double u[3], double v[3], double w[3], double outputBounds[6]);
119
121
125 void SetMinPoint(double x, double y, double z);
126 void SetMinPoint(double p[3]);
128
130
134 void SetMaxPoint(double x, double y, double z);
135 void SetMaxPoint(double p[3]);
137
139
143 int IsValid() const;
144 static int IsValid(const double bounds[6]);
146
148
152 void AddPoint(double p[3]);
153 void AddPoint(double px, double py, double pz);
155
160 void AddBox(const vtkBoundingBox& bbox);
161
166 void AddBounds(const double bounds[]);
167
171 bool IsSubsetOf(const vtkBoundingBox& bbox) const;
172
178 int IntersectBox(const vtkBoundingBox& bbox);
179
183 int Intersects(const vtkBoundingBox& bbox) const;
184
190 bool IntersectPlane(double origin[3], double normal[3]);
191
196 bool IntersectsSphere(double center[3], double squaredRadius) const;
197
202 bool IntersectsLine(const double p1[3], const double p2[3]) const;
203
208
213 int Contains(const vtkBoundingBox& bbox) const;
214
216
219 void GetBounds(double bounds[6]) const;
220 void GetBounds(
221 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const;
223
227 double GetBound(int i) const;
228
230
233 const double* GetMinPoint() const VTK_SIZEHINT(3);
234 void GetMinPoint(double& x, double& y, double& z) const;
235 void GetMinPoint(double x[3]) const;
237
239
242 const double* GetMaxPoint() const VTK_SIZEHINT(3);
243 void GetMaxPoint(double& x, double& y, double& z) const;
244 void GetMaxPoint(double x[3]) const;
246
251 void GetCorner(int corner, double p[3]) const;
252
254
257 vtkTypeBool ContainsPoint(const double p[3]) const;
258 vtkTypeBool ContainsPoint(double px, double py, double pz) const;
259 template <class PointT>
260 bool ContainsPoint(const PointT& p) const;
262
266 void GetCenter(double center[3]) const;
267
271 void GetLengths(double lengths[3]) const;
272
276 double GetLength(int i) const;
277
281 double GetMaxLength() const;
282
287 double GetDiagonalLength() const;
288
290
298 void Inflate(double delta);
299 void Inflate(double deltaX, double deltaY, double deltaZ);
300 void Inflate();
302
304
310 void Scale(double s[3]);
311 void Scale(double sx, double sy, double sz);
313
315
320 void ScaleAboutCenter(double s);
321 void ScaleAboutCenter(double s[3]);
322 void ScaleAboutCenter(double sx, double sy, double sz);
324
335 vtkIdType ComputeDivisions(vtkIdType totalBins, double bounds[6], int divs[3]) const;
336
341 static void ClampDivisions(vtkIdType targetBins, int divs[3]);
342
346 void Reset();
347
348protected:
349 double MinPnt[3], MaxPnt[3];
350};
351
352inline void vtkBoundingBox::Reset()
353{
354 this->MinPnt[0] = this->MinPnt[1] = this->MinPnt[2] = VTK_DOUBLE_MAX;
355 this->MaxPnt[0] = this->MaxPnt[1] = this->MaxPnt[2] = VTK_DOUBLE_MIN;
356}
357
359 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const
360{
361 xMin = this->MinPnt[0];
362 xMax = this->MaxPnt[0];
363 yMin = this->MinPnt[1];
364 yMax = this->MaxPnt[1];
365 zMin = this->MinPnt[2];
366 zMax = this->MaxPnt[2];
367}
368
369inline double vtkBoundingBox::GetBound(int i) const
370{
371 // If i is odd then when are returning a part of the max bounds
372 // else part of the min bounds is requested. The exact component
373 // needed is i /2 (or i right shifted by 1
374 return ((i & 0x1) ? this->MaxPnt[i >> 1] : this->MinPnt[i >> 1]);
375}
376
377inline const double* vtkBoundingBox::GetMinPoint() const
378{
379 return this->MinPnt;
380}
381
382inline void vtkBoundingBox::GetMinPoint(double x[3]) const
383{
384 x[0] = this->MinPnt[0];
385 x[1] = this->MinPnt[1];
386 x[2] = this->MinPnt[2];
387}
388
389inline const double* vtkBoundingBox::GetMaxPoint() const
390{
391 return this->MaxPnt;
392}
393
394inline void vtkBoundingBox::GetMaxPoint(double x[3]) const
395{
396 x[0] = this->MaxPnt[0];
397 x[1] = this->MaxPnt[1];
398 x[2] = this->MaxPnt[2];
399}
400
401inline int vtkBoundingBox::IsValid() const
402{
403 return ((this->MinPnt[0] <= this->MaxPnt[0]) && (this->MinPnt[1] <= this->MaxPnt[1]) &&
404 (this->MinPnt[2] <= this->MaxPnt[2]));
405}
406
407inline int vtkBoundingBox::IsValid(const double bounds[6])
408{
409 return (bounds[0] <= bounds[1] && bounds[2] <= bounds[3] && bounds[4] <= bounds[5]);
410}
411
412inline double vtkBoundingBox::GetLength(int i) const
413{
414 return this->MaxPnt[i] - this->MinPnt[i];
415}
416
417inline void vtkBoundingBox::GetLengths(double lengths[3]) const
418{
419 lengths[0] = this->GetLength(0);
420 lengths[1] = this->GetLength(1);
421 lengths[2] = this->GetLength(2);
422}
423
424inline void vtkBoundingBox::GetCenter(double center[3]) const
425{
426 center[0] = 0.5 * (this->MaxPnt[0] + this->MinPnt[0]);
427 center[1] = 0.5 * (this->MaxPnt[1] + this->MinPnt[1]);
428 center[2] = 0.5 * (this->MaxPnt[2] + this->MinPnt[2]);
429}
430
431inline bool vtkBoundingBox::IsSubsetOf(const vtkBoundingBox& bbox) const
432{
433 const double* bboxMaxPnt = bbox.GetMaxPoint();
434 const double* bboxMinPnt = bbox.GetMinPoint();
435 return this->MaxPnt[0] < bboxMaxPnt[0] && this->MinPnt[0] > bboxMinPnt[0] &&
436 this->MaxPnt[1] < bboxMaxPnt[1] && this->MinPnt[1] > bboxMinPnt[1] &&
437 this->MaxPnt[2] < bboxMaxPnt[2] && this->MinPnt[2] > bboxMinPnt[2];
438}
439
440inline void vtkBoundingBox::SetBounds(const double bounds[6])
441{
442 this->SetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
443}
444
445inline void vtkBoundingBox::GetBounds(double bounds[6]) const
446{
447 this->GetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
448}
449
451{
452 this->Reset();
453}
454
455inline vtkBoundingBox::vtkBoundingBox(const double bounds[6])
456{
457 this->Reset();
458 this->SetBounds(bounds);
459}
460
462 double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
463{
464 this->Reset();
465 this->SetBounds(xMin, xMax, yMin, yMax, zMin, zMax);
466}
467
469{
470 this->MinPnt[0] = bbox.MinPnt[0];
471 this->MinPnt[1] = bbox.MinPnt[1];
472 this->MinPnt[2] = bbox.MinPnt[2];
473
474 this->MaxPnt[0] = bbox.MaxPnt[0];
475 this->MaxPnt[1] = bbox.MaxPnt[1];
476 this->MaxPnt[2] = bbox.MaxPnt[2];
477}
478
480{
481 this->MinPnt[0] = bbox.MinPnt[0];
482 this->MinPnt[1] = bbox.MinPnt[1];
483 this->MinPnt[2] = bbox.MinPnt[2];
484
485 this->MaxPnt[0] = bbox.MaxPnt[0];
486 this->MaxPnt[1] = bbox.MaxPnt[1];
487 this->MaxPnt[2] = bbox.MaxPnt[2];
488 return *this;
489}
490
491inline bool vtkBoundingBox::operator==(const vtkBoundingBox& bbox) const
492{
493 return ((this->MinPnt[0] == bbox.MinPnt[0]) && (this->MinPnt[1] == bbox.MinPnt[1]) &&
494 (this->MinPnt[2] == bbox.MinPnt[2]) && (this->MaxPnt[0] == bbox.MaxPnt[0]) &&
495 (this->MaxPnt[1] == bbox.MaxPnt[1]) && (this->MaxPnt[2] == bbox.MaxPnt[2]));
496}
497
498inline bool vtkBoundingBox::operator!=(const vtkBoundingBox& bbox) const
499{
500 return !((*this) == bbox);
501}
502
503inline void vtkBoundingBox::SetMinPoint(double p[3])
504{
505 this->SetMinPoint(p[0], p[1], p[2]);
506}
507
508inline void vtkBoundingBox::SetMaxPoint(double p[3])
509{
510 this->SetMaxPoint(p[0], p[1], p[2]);
511}
512
513inline void vtkBoundingBox::GetMinPoint(double& x, double& y, double& z) const
514{
515 x = this->MinPnt[0];
516 y = this->MinPnt[1];
517 z = this->MinPnt[2];
518}
519
520inline void vtkBoundingBox::GetMaxPoint(double& x, double& y, double& z) const
521{
522 x = this->MaxPnt[0];
523 y = this->MaxPnt[1];
524 z = this->MaxPnt[2];
525}
526
527inline vtkTypeBool vtkBoundingBox::ContainsPoint(double px, double py, double pz) const
528{
529 if ((px < this->MinPnt[0]) || (px > this->MaxPnt[0]))
530 {
531 return 0;
532 }
533 if ((py < this->MinPnt[1]) || (py > this->MaxPnt[1]))
534 {
535 return 0;
536 }
537 if ((pz < this->MinPnt[2]) || (pz > this->MaxPnt[2]))
538 {
539 return 0;
540 }
541 return 1;
542}
543
544inline vtkTypeBool vtkBoundingBox::ContainsPoint(const double p[3]) const
545{
546 return this->ContainsPoint(p[0], p[1], p[2]);
547}
548
549template <class PointT>
550inline bool vtkBoundingBox::ContainsPoint(const PointT& p) const
551{
552 return this->ContainsPoint(p[0], p[1], p[2]);
553}
554
555inline void vtkBoundingBox::GetCorner(int corner, double p[3]) const
556{
557 if ((corner < 0) || (corner > 7))
558 {
559 p[0] = VTK_DOUBLE_MAX;
560 p[1] = VTK_DOUBLE_MAX;
561 p[2] = VTK_DOUBLE_MAX;
562 return; // out of bounds
563 }
564
565 int ix = (corner & 1) ? 1 : 0; // 0,1,0,1,0,1,0,1
566 int iy = ((corner >> 1) & 1) ? 1 : 0; // 0,0,1,1,0,0,1,1
567 int iz = (corner >> 2) ? 1 : 0; // 0,0,0,0,1,1,1,1
568
569 const double* pts[2] = { this->MinPnt, this->MaxPnt };
570 p[0] = pts[ix][0];
571 p[1] = pts[iy][1];
572 p[2] = pts[iz][2];
573}
574
575#endif
576// VTK-HeaderTest-Exclude: vtkBoundingBox.h
Fast, simple class for representing and operating on 3D bounds.
int IntersectBox(const vtkBoundingBox &bbox)
Intersect this box with bbox.
const double * GetMinPoint() const
Get the minimum point of the bounding box.
void SetBounds(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
Set the bounds explicitly of the box (using the VTK convention for representing a bounding box).
void AddBox(const vtkBoundingBox &bbox)
Change the bounding box to be the union of itself and the specified bbox.
void AddBounds(const double bounds[])
Adjust the bounding box so it contains the specified bounds (defined by the VTK representation (xmin,...
int Contains(const vtkBoundingBox &bbox) const
Returns 1 if the min and max points of bbox are contained within the bounds of the specified box,...
int IsValid() const
Returns 1 if the bounds have been set and 0 if the box is in its initialized state which is an invert...
int Intersects(const vtkBoundingBox &bbox) const
Returns 1 if the boxes intersect else returns 0.
bool operator!=(const vtkBoundingBox &bbox) const
Equality operator.
void AddPoint(double px, double py, double pz)
Change bounding box so it includes the point p.
int ComputeInnerDimension() const
Returns the inner dimension of the bounding box.
void GetCorner(int corner, double p[3]) const
Get the ith corner of the bounding box.
void ComputeBounds(vtkPoints *pts)
Compute the bounding box from an array of vtkPoints.
bool IsSubsetOf(const vtkBoundingBox &bbox) const
Returns true if this instance is entirely contained by bbox.
static void ComputeBounds(vtkPoints *pts, double bounds[6])
Compute the bounding box from an array of vtkPoints.
bool IntersectsSphere(double center[3], double squaredRadius) const
Intersect this box with a sphere.
void SetMaxPoint(double x, double y, double z)
Set the maximum point of the bounding box - if the max point is less than the min point then the min ...
bool IntersectPlane(double origin[3], double normal[3])
Intersect this box with the half space defined by plane.
bool IntersectsLine(const double p1[3], const double p2[3]) const
Returns true if any part of segment [p1,p2] lies inside the bounding box, as well as on its boundarie...
static void ComputeLocalBounds(vtkPoints *points, double u[3], double v[3], double w[3], double outputBounds[6])
Compute local bounds.
void GetCenter(double center[3]) const
Get the center of the bounding box.
void AddPoint(double p[3])
Change bounding box so it includes the point p.
double GetLength(int i) const
Return the length of the bounding box in the ith direction.
bool operator==(const vtkBoundingBox &bbox) const
Equality operator.
vtkTypeBool ContainsPoint(const double p[3]) const
Returns 1 if the point is contained in the box else 0.
vtkBoundingBox()
Construct a bounding box with the min point set to VTK_DOUBLE_MAX and the max point set to VTK_DOUBLE...
double MinPnt[3]
double MaxPnt[3]
void GetLengths(double lengths[3]) const
Get the length of each side of the box.
void ComputeBounds(vtkPoints *pts, unsigned char *ptUses)
Compute the bounding box from an array of vtkPoints.
static void ComputeBounds(vtkPoints *pts, const unsigned char *ptUses, double bounds[6])
Compute the bounding box from an array of vtkPoints.
void SetBounds(const double bounds[6])
Set the bounds explicitly of the box (using the VTK convention for representing a bounding box).
const double * GetMaxPoint() const
Get the maximum point of the bounding box.
double GetBound(int i) const
Return the ith bounds of the box (defined by VTK style).
static void ComputeBounds(vtkPoints *pts, const std::atomic< unsigned char > *ptUses, double bounds[6])
Compute the bounding box from an array of vtkPoints.
void GetBounds(double bounds[6]) const
Get the bounds of the box (defined by VTK style).
void SetMinPoint(double x, double y, double z)
Set the minimum point of the bounding box - if the min point is greater than the max point then the m...
vtkBoundingBox & operator=(const vtkBoundingBox &bbox)
Assignment Operator.
represent and manipulate 3D points
Definition: vtkPoints.h:43
void GetBounds(T a, double bds[6])
@ points
Definition: vtkX3D.h:452
@ center
Definition: vtkX3D.h:236
int vtkTypeBool
Definition: vtkABI.h:69
int vtkIdType
Definition: vtkType.h:332
#define VTK_DOUBLE_MIN
Definition: vtkType.h:164
#define VTK_DOUBLE_MAX
Definition: vtkType.h:165
VTKCOMMONCORE_EXPORT bool operator!=(const vtkUnicodeString &lhs, const vtkUnicodeString &rhs)
VTKCOMMONCORE_EXPORT bool operator==(const vtkUnicodeString &lhs, const vtkUnicodeString &rhs)
#define VTK_SIZEHINT(...)