VTK  9.5.20251018
vtkBoundingBox.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
48#ifndef vtkBoundingBox_h
49#define vtkBoundingBox_h
50#include "vtkCommonDataModelModule.h" // For export macro
51#include "vtkSystemIncludes.h"
52#include <atomic> // For threaded bounding box computation
53
54VTK_ABI_NAMESPACE_BEGIN
55class vtkPoints;
56
57class VTKCOMMONDATAMODEL_EXPORT vtkBoundingBox
58{
59public:
61
69 vtkBoundingBox(const double bounds[6]);
73 vtkBoundingBox(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
77 vtkBoundingBox(double center[3], double delta);
79
83 vtkBoundingBox(const vtkBoundingBox& bbox);
84
88 vtkBoundingBox& operator=(const vtkBoundingBox& bbox);
89
91
94 bool operator==(const vtkBoundingBox& bbox) const;
95 bool operator!=(const vtkBoundingBox& bbox) const;
97
99
103 void SetBounds(const double bounds[6]);
104 void SetBounds(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
106
108
115 static void ComputeBounds(vtkPoints* pts, double bounds[6]);
116 static void ComputeBounds(vtkPoints* pts, const unsigned char* ptUses, double bounds[6]);
117 static void ComputeBounds(
118 vtkPoints* pts, const std::atomic<unsigned char>* ptUses, double bounds[6]);
119 template <typename TIter>
120 static void ComputeBounds(vtkPoints* pts, TIter ptIds, vtkIdType numPointIds, double bounds[6]);
122 {
123 double bds[6];
125 this->MinPnt[0] = bds[0];
126 this->MinPnt[1] = bds[2];
127 this->MinPnt[2] = bds[4];
128 this->MaxPnt[0] = bds[1];
129 this->MaxPnt[1] = bds[3];
130 this->MaxPnt[2] = bds[5];
131 }
132 void ComputeBounds(vtkPoints* pts, unsigned char* ptUses)
133 {
134 double bds[6];
135 vtkBoundingBox::ComputeBounds(pts, ptUses, bds);
136 this->MinPnt[0] = bds[0];
137 this->MinPnt[1] = bds[2];
138 this->MinPnt[2] = bds[4];
139 this->MaxPnt[0] = bds[1];
140 this->MaxPnt[1] = bds[3];
141 this->MaxPnt[2] = bds[5];
142 }
144
146
151 vtkPoints* points, double u[3], double v[3], double w[3], double outputBounds[6]);
153
155
159 void SetMinPoint(double x, double y, double z);
160 void SetMinPoint(double p[3]);
162
164
168 void SetMaxPoint(double x, double y, double z);
169 void SetMaxPoint(double p[3]);
171
173
177 int IsValid() const;
178 static int IsValid(const double bounds[6]);
180
182
186 void AddPoint(double p[3]);
187 void AddPoint(double px, double py, double pz);
189
194 void AddBox(const vtkBoundingBox& bbox);
195
200 void AddBounds(const double bounds[6]);
201
205 bool IsSubsetOf(const vtkBoundingBox& bbox) const;
206
212 int IntersectBox(const vtkBoundingBox& bbox);
213
217 int Intersects(const vtkBoundingBox& bbox) const;
218
224 bool IntersectPlane(double origin[3], double normal[3]);
225
230 bool IntersectsSphere(double center[3], double squaredRadius) const;
231
236 bool IntersectsLine(const double p1[3], const double p2[3]) const;
237
242
247 int Contains(const vtkBoundingBox& bbox) const;
248
265 static bool ContainsLine(const double x[3], const double s[3], const double lineEnd[3], double& t,
266 double xInt[3], int& plane);
267
269
272 void GetBounds(double bounds[6]) const;
273 void GetBounds(
274 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const;
276
280 double GetBound(int i) const;
281
283
286 const double* GetMinPoint() const VTK_SIZEHINT(3);
287 void GetMinPoint(double& x, double& y, double& z) const;
288 void GetMinPoint(double x[3]) const;
290
292
295 const double* GetMaxPoint() const VTK_SIZEHINT(3);
296 void GetMaxPoint(double& x, double& y, double& z) const;
297 void GetMaxPoint(double x[3]) const;
299
304 void GetCorner(int corner, double p[3]) const;
305
307
310 vtkTypeBool ContainsPoint(const double p[3]) const;
311 vtkTypeBool ContainsPoint(double px, double py, double pz) const;
312 template <class PointT>
313 bool ContainsPoint(const PointT& p) const;
315
319 void GetCenter(double center[3]) const;
320
324 void GetLengths(double lengths[3]) const;
325
329 double GetLength(int i) const;
330
334 double GetMaxLength() const;
335
337
341 double GetDiagonalLength2() const;
342 double GetDiagonalLength() const;
344
346
357 void Inflate(double delta);
358 void Inflate(double deltaX, double deltaY, double deltaZ);
359 void Inflate();
360 void InflateSlice(double delta);
362
364
370 void Scale(double s[3]);
371 void Scale(double sx, double sy, double sz);
373
375
380 void ScaleAboutCenter(double s);
381 void ScaleAboutCenter(double s[3]);
382 void ScaleAboutCenter(double sx, double sy, double sz);
384
395 vtkIdType ComputeDivisions(vtkIdType totalBins, double bounds[6], int divs[3]) const;
396
401 static void ClampDivisions(vtkIdType targetBins, int divs[3]);
402
406 void Reset();
407
412 void ClampPoint(double point[3]);
413
420 void GetDistance(double point[3], double distance[3]);
421
426 void Translate(double motion[3]);
427
428protected:
429 double MinPnt[3], MaxPnt[3];
430};
431
432inline void vtkBoundingBox::Reset()
433{
434 this->MinPnt[0] = this->MinPnt[1] = this->MinPnt[2] = VTK_DOUBLE_MAX;
435 this->MaxPnt[0] = this->MaxPnt[1] = this->MaxPnt[2] = VTK_DOUBLE_MIN;
436}
437
439 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const
440{
441 xMin = this->MinPnt[0];
442 xMax = this->MaxPnt[0];
443 yMin = this->MinPnt[1];
444 yMax = this->MaxPnt[1];
445 zMin = this->MinPnt[2];
446 zMax = this->MaxPnt[2];
447}
448
449inline double vtkBoundingBox::GetBound(int i) const
450{
451 // If i is odd then when are returning a part of the max bounds
452 // else part of the min bounds is requested. The exact component
453 // needed is i /2 (or i right shifted by 1
454 return ((i & 0x1) ? this->MaxPnt[i >> 1] : this->MinPnt[i >> 1]);
455}
456
457inline const double* vtkBoundingBox::GetMinPoint() const
458{
459 return this->MinPnt;
460}
461
462inline void vtkBoundingBox::GetMinPoint(double x[3]) const
463{
464 x[0] = this->MinPnt[0];
465 x[1] = this->MinPnt[1];
466 x[2] = this->MinPnt[2];
467}
468
469inline const double* vtkBoundingBox::GetMaxPoint() const
470{
471 return this->MaxPnt;
472}
473
474inline void vtkBoundingBox::GetMaxPoint(double x[3]) const
475{
476 x[0] = this->MaxPnt[0];
477 x[1] = this->MaxPnt[1];
478 x[2] = this->MaxPnt[2];
479}
480
481inline int vtkBoundingBox::IsValid() const
482{
483 return ((this->MinPnt[0] <= this->MaxPnt[0]) && (this->MinPnt[1] <= this->MaxPnt[1]) &&
484 (this->MinPnt[2] <= this->MaxPnt[2]));
485}
486
487inline int vtkBoundingBox::IsValid(const double bounds[6])
488{
489 return (bounds[0] <= bounds[1] && bounds[2] <= bounds[3] && bounds[4] <= bounds[5]);
490}
491
492inline double vtkBoundingBox::GetLength(int i) const
493{
494 return this->MaxPnt[i] - this->MinPnt[i];
495}
496
497inline void vtkBoundingBox::GetLengths(double lengths[3]) const
498{
499 lengths[0] = this->GetLength(0);
500 lengths[1] = this->GetLength(1);
501 lengths[2] = this->GetLength(2);
502}
503
504inline void vtkBoundingBox::GetCenter(double center[3]) const
505{
506 center[0] = 0.5 * (this->MaxPnt[0] + this->MinPnt[0]);
507 center[1] = 0.5 * (this->MaxPnt[1] + this->MinPnt[1]);
508 center[2] = 0.5 * (this->MaxPnt[2] + this->MinPnt[2]);
509}
510
511inline bool vtkBoundingBox::IsSubsetOf(const vtkBoundingBox& bbox) const
512{
513 const double* bboxMaxPnt = bbox.GetMaxPoint();
514 const double* bboxMinPnt = bbox.GetMinPoint();
515 return this->MaxPnt[0] < bboxMaxPnt[0] && this->MinPnt[0] > bboxMinPnt[0] &&
516 this->MaxPnt[1] < bboxMaxPnt[1] && this->MinPnt[1] > bboxMinPnt[1] &&
517 this->MaxPnt[2] < bboxMaxPnt[2] && this->MinPnt[2] > bboxMinPnt[2];
518}
519
520inline void vtkBoundingBox::SetBounds(const double bounds[6])
521{
522 this->SetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
523}
524
525inline void vtkBoundingBox::GetBounds(double bounds[6]) const
526{
527 this->GetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
528}
529
531{
532 this->Reset();
533}
534
535inline vtkBoundingBox::vtkBoundingBox(const double bounds[6])
536{
537 this->Reset();
538 this->SetBounds(bounds);
539}
540
542 double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
543{
544 this->Reset();
545 this->SetBounds(xMin, xMax, yMin, yMax, zMin, zMax);
546}
547
549{
550 this->MinPnt[0] = bbox.MinPnt[0];
551 this->MinPnt[1] = bbox.MinPnt[1];
552 this->MinPnt[2] = bbox.MinPnt[2];
553
554 this->MaxPnt[0] = bbox.MaxPnt[0];
555 this->MaxPnt[1] = bbox.MaxPnt[1];
556 this->MaxPnt[2] = bbox.MaxPnt[2];
557}
558
559inline vtkBoundingBox::vtkBoundingBox(double center[3], double delta)
560{
561 this->Reset();
562 this->AddPoint(center);
563 this->Inflate(delta);
564}
565
567{
568 this->MinPnt[0] = bbox.MinPnt[0];
569 this->MinPnt[1] = bbox.MinPnt[1];
570 this->MinPnt[2] = bbox.MinPnt[2];
571
572 this->MaxPnt[0] = bbox.MaxPnt[0];
573 this->MaxPnt[1] = bbox.MaxPnt[1];
574 this->MaxPnt[2] = bbox.MaxPnt[2];
575 return *this;
576}
577
578inline bool vtkBoundingBox::operator==(const vtkBoundingBox& bbox) const
579{
580 return ((this->MinPnt[0] == bbox.MinPnt[0]) && (this->MinPnt[1] == bbox.MinPnt[1]) &&
581 (this->MinPnt[2] == bbox.MinPnt[2]) && (this->MaxPnt[0] == bbox.MaxPnt[0]) &&
582 (this->MaxPnt[1] == bbox.MaxPnt[1]) && (this->MaxPnt[2] == bbox.MaxPnt[2]));
583}
584
585inline bool vtkBoundingBox::operator!=(const vtkBoundingBox& bbox) const
586{
587 return !((*this) == bbox);
588}
589
590inline void vtkBoundingBox::SetMinPoint(double p[3])
591{
592 this->SetMinPoint(p[0], p[1], p[2]);
593}
594
595inline void vtkBoundingBox::SetMaxPoint(double p[3])
596{
597 this->SetMaxPoint(p[0], p[1], p[2]);
598}
599
600inline void vtkBoundingBox::GetMinPoint(double& x, double& y, double& z) const
601{
602 x = this->MinPnt[0];
603 y = this->MinPnt[1];
604 z = this->MinPnt[2];
605}
606
607inline void vtkBoundingBox::GetMaxPoint(double& x, double& y, double& z) const
608{
609 x = this->MaxPnt[0];
610 y = this->MaxPnt[1];
611 z = this->MaxPnt[2];
612}
613
614inline vtkTypeBool vtkBoundingBox::ContainsPoint(double px, double py, double pz) const
615{
616 if ((px < this->MinPnt[0]) || (px > this->MaxPnt[0]))
617 {
618 return 0;
619 }
620 if ((py < this->MinPnt[1]) || (py > this->MaxPnt[1]))
621 {
622 return 0;
623 }
624 if ((pz < this->MinPnt[2]) || (pz > this->MaxPnt[2]))
625 {
626 return 0;
627 }
628 return 1;
629}
630
631inline vtkTypeBool vtkBoundingBox::ContainsPoint(const double p[3]) const
632{
633 return this->ContainsPoint(p[0], p[1], p[2]);
634}
635
636template <class PointT>
637inline bool vtkBoundingBox::ContainsPoint(const PointT& p) const
638{
639 return this->ContainsPoint(p[0], p[1], p[2]);
640}
641
642inline void vtkBoundingBox::GetCorner(int corner, double p[3]) const
643{
644 if ((corner < 0) || (corner > 7))
645 {
646 p[0] = VTK_DOUBLE_MAX;
647 p[1] = VTK_DOUBLE_MAX;
648 p[2] = VTK_DOUBLE_MAX;
649 return; // out of bounds
650 }
651
652 int ix = (corner & 1) ? 1 : 0; // 0,1,0,1,0,1,0,1
653 int iy = ((corner >> 1) & 1) ? 1 : 0; // 0,0,1,1,0,0,1,1
654 int iz = (corner >> 2) ? 1 : 0; // 0,0,0,0,1,1,1,1
655
656 const double* pts[2] = { this->MinPnt, this->MaxPnt };
657 p[0] = pts[ix][0];
658 p[1] = pts[iy][1];
659 p[2] = pts[iz][2];
660}
661
662VTK_ABI_NAMESPACE_END
663#endif
664// VTK-HeaderTest-Exclude: vtkBoundingBox.h
Fast, simple class for representing and operating on 3D bounds.
static bool ContainsLine(const double x[3], const double s[3], const double lineEnd[3], double &t, double xInt[3], int &plane)
A specialized, performant method to compute the containment of a finite line emanating from the cente...
void AddBounds(const double bounds[6])
Adjust the bounding box so it contains the specified bounds (defined by the VTK representation (xmin,...
static void ComputeBounds(vtkPoints *pts, TIter ptIds, vtkIdType numPointIds, double bounds[6])
Compute the bounding box from an array of vtkPoints.
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.
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...
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:139
int vtkTypeBool
Definition vtkABI.h:64
bool VTKCOMMONCORE_EXPORT operator==(const std::string &a, const vtkStringToken &b)
bool VTKCOMMONCORE_EXPORT operator!=(const std::string &a, const vtkStringToken &b)
int vtkIdType
Definition vtkType.h:367
#define VTK_DOUBLE_MIN
Definition vtkType.h:205
#define VTK_DOUBLE_MAX
Definition vtkType.h:206
#define VTK_SIZEHINT(...)