VTK  9.5.2
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
18
19#ifndef vtkBoundingBox_h
20#define vtkBoundingBox_h
21#include "vtkCommonDataModelModule.h" // For export macro
22#include "vtkSystemIncludes.h"
23#include <atomic> // For threaded bounding box computation
24
25VTK_ABI_NAMESPACE_BEGIN
26class vtkPoints;
27
28class VTKCOMMONDATAMODEL_EXPORT vtkBoundingBox
29{
30public:
32
40 vtkBoundingBox(const double bounds[6]);
44 vtkBoundingBox(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
48 vtkBoundingBox(double center[3], double delta);
50
54 vtkBoundingBox(const vtkBoundingBox& bbox);
55
60
62
65 bool operator==(const vtkBoundingBox& bbox) const;
66 bool operator!=(const vtkBoundingBox& bbox) const;
68
70
74 void SetBounds(const double bounds[6]);
75 void SetBounds(double xMin, double xMax, double yMin, double yMax, double zMin, double zMax);
77
79
86 static void ComputeBounds(vtkPoints* pts, double bounds[6]);
87 static void ComputeBounds(vtkPoints* pts, const unsigned char* ptUses, double bounds[6]);
88 static void ComputeBounds(
89 vtkPoints* pts, const std::atomic<unsigned char>* ptUses, double bounds[6]);
90 static void ComputeBounds(
91 vtkPoints* pts, const long long* ptIds, long long numPointIds, double bounds[6]);
92 static void ComputeBounds(vtkPoints* pts, const long* ptIds, long numPointIds, double bounds[6]);
93 static void ComputeBounds(vtkPoints* pts, const int* ptIds, int numPointIds, double bounds[6]);
95 {
96 double bds[6];
98 this->MinPnt[0] = bds[0];
99 this->MinPnt[1] = bds[2];
100 this->MinPnt[2] = bds[4];
101 this->MaxPnt[0] = bds[1];
102 this->MaxPnt[1] = bds[3];
103 this->MaxPnt[2] = bds[5];
104 }
105 void ComputeBounds(vtkPoints* pts, unsigned char* ptUses)
106 {
107 double bds[6];
108 vtkBoundingBox::ComputeBounds(pts, ptUses, bds);
109 this->MinPnt[0] = bds[0];
110 this->MinPnt[1] = bds[2];
111 this->MinPnt[2] = bds[4];
112 this->MaxPnt[0] = bds[1];
113 this->MaxPnt[1] = bds[3];
114 this->MaxPnt[2] = bds[5];
115 }
116
117
119
124 vtkPoints* points, double u[3], double v[3], double w[3], double outputBounds[6]);
126
128
132 void SetMinPoint(double x, double y, double z);
133 void SetMinPoint(double p[3]);
135
137
141 void SetMaxPoint(double x, double y, double z);
142 void SetMaxPoint(double p[3]);
144
146
150 int IsValid() const;
151 static int IsValid(const double bounds[6]);
153
155
159 void AddPoint(double p[3]);
160 void AddPoint(double px, double py, double pz);
162
167 void AddBox(const vtkBoundingBox& bbox);
168
173 void AddBounds(const double bounds[6]);
174
178 bool IsSubsetOf(const vtkBoundingBox& bbox) const;
179
185 int IntersectBox(const vtkBoundingBox& bbox);
186
190 int Intersects(const vtkBoundingBox& bbox) const;
191
197 bool IntersectPlane(double origin[3], double normal[3]);
198
203 bool IntersectsSphere(double center[3], double squaredRadius) const;
204
209 bool IntersectsLine(const double p1[3], const double p2[3]) const;
210
215
220 int Contains(const vtkBoundingBox& bbox) const;
221
238 static bool ContainsLine(const double x[3], const double s[3], const double lineEnd[3], double& t,
239 double xInt[3], int& plane);
240
242
245 void GetBounds(double bounds[6]) const;
246 void GetBounds(
247 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const;
249
253 double GetBound(int i) const;
254
256
259 const double* GetMinPoint() const VTK_SIZEHINT(3);
260 void GetMinPoint(double& x, double& y, double& z) const;
261 void GetMinPoint(double x[3]) const;
263
265
268 const double* GetMaxPoint() const VTK_SIZEHINT(3);
269 void GetMaxPoint(double& x, double& y, double& z) const;
270 void GetMaxPoint(double x[3]) const;
272
277 void GetCorner(int corner, double p[3]) const;
278
280
283 vtkTypeBool ContainsPoint(const double p[3]) const;
284 vtkTypeBool ContainsPoint(double px, double py, double pz) const;
285 template <class PointT>
286 bool ContainsPoint(const PointT& p) const;
288
292 void GetCenter(double center[3]) const;
293
297 void GetLengths(double lengths[3]) const;
298
302 double GetLength(int i) const;
303
307 double GetMaxLength() const;
308
310
314 double GetDiagonalLength2() const;
315 double GetDiagonalLength() const;
317
319
330 void Inflate(double delta);
331 void Inflate(double deltaX, double deltaY, double deltaZ);
332 void Inflate();
333 void InflateSlice(double delta);
335
337
343 void Scale(double s[3]);
344 void Scale(double sx, double sy, double sz);
346
348
353 void ScaleAboutCenter(double s);
354 void ScaleAboutCenter(double s[3]);
355 void ScaleAboutCenter(double sx, double sy, double sz);
357
368 vtkIdType ComputeDivisions(vtkIdType totalBins, double bounds[6], int divs[3]) const;
369
374 static void ClampDivisions(vtkIdType targetBins, int divs[3]);
375
379 void Reset();
380
385 void ClampPoint(double point[3]);
386
393 void GetDistance(double point[3], double distance[3]);
394
399 void Translate(double motion[3]);
400
401protected:
402 double MinPnt[3], MaxPnt[3];
403};
404
405inline void vtkBoundingBox::Reset()
406{
407 this->MinPnt[0] = this->MinPnt[1] = this->MinPnt[2] = VTK_DOUBLE_MAX;
408 this->MaxPnt[0] = this->MaxPnt[1] = this->MaxPnt[2] = VTK_DOUBLE_MIN;
409}
410
412 double& xMin, double& xMax, double& yMin, double& yMax, double& zMin, double& zMax) const
413{
414 xMin = this->MinPnt[0];
415 xMax = this->MaxPnt[0];
416 yMin = this->MinPnt[1];
417 yMax = this->MaxPnt[1];
418 zMin = this->MinPnt[2];
419 zMax = this->MaxPnt[2];
420}
421
422inline double vtkBoundingBox::GetBound(int i) const
423{
424 // If i is odd then when are returning a part of the max bounds
425 // else part of the min bounds is requested. The exact component
426 // needed is i /2 (or i right shifted by 1
427 return ((i & 0x1) ? this->MaxPnt[i >> 1] : this->MinPnt[i >> 1]);
428}
429
430inline const double* vtkBoundingBox::GetMinPoint() const
431{
432 return this->MinPnt;
433}
434
435inline void vtkBoundingBox::GetMinPoint(double x[3]) const
436{
437 x[0] = this->MinPnt[0];
438 x[1] = this->MinPnt[1];
439 x[2] = this->MinPnt[2];
440}
441
442inline const double* vtkBoundingBox::GetMaxPoint() const
443{
444 return this->MaxPnt;
445}
446
447inline void vtkBoundingBox::GetMaxPoint(double x[3]) const
448{
449 x[0] = this->MaxPnt[0];
450 x[1] = this->MaxPnt[1];
451 x[2] = this->MaxPnt[2];
452}
453
454inline int vtkBoundingBox::IsValid() const
455{
456 return ((this->MinPnt[0] <= this->MaxPnt[0]) && (this->MinPnt[1] <= this->MaxPnt[1]) &&
457 (this->MinPnt[2] <= this->MaxPnt[2]));
458}
459
460inline int vtkBoundingBox::IsValid(const double bounds[6])
461{
462 return (bounds[0] <= bounds[1] && bounds[2] <= bounds[3] && bounds[4] <= bounds[5]);
463}
464
465inline double vtkBoundingBox::GetLength(int i) const
466{
467 return this->MaxPnt[i] - this->MinPnt[i];
468}
469
470inline void vtkBoundingBox::GetLengths(double lengths[3]) const
471{
472 lengths[0] = this->GetLength(0);
473 lengths[1] = this->GetLength(1);
474 lengths[2] = this->GetLength(2);
475}
476
477inline void vtkBoundingBox::GetCenter(double center[3]) const
478{
479 center[0] = 0.5 * (this->MaxPnt[0] + this->MinPnt[0]);
480 center[1] = 0.5 * (this->MaxPnt[1] + this->MinPnt[1]);
481 center[2] = 0.5 * (this->MaxPnt[2] + this->MinPnt[2]);
482}
483
484inline bool vtkBoundingBox::IsSubsetOf(const vtkBoundingBox& bbox) const
485{
486 const double* bboxMaxPnt = bbox.GetMaxPoint();
487 const double* bboxMinPnt = bbox.GetMinPoint();
488 return this->MaxPnt[0] < bboxMaxPnt[0] && this->MinPnt[0] > bboxMinPnt[0] &&
489 this->MaxPnt[1] < bboxMaxPnt[1] && this->MinPnt[1] > bboxMinPnt[1] &&
490 this->MaxPnt[2] < bboxMaxPnt[2] && this->MinPnt[2] > bboxMinPnt[2];
491}
492
493inline void vtkBoundingBox::SetBounds(const double bounds[6])
494{
495 this->SetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
496}
497
498inline void vtkBoundingBox::GetBounds(double bounds[6]) const
499{
500 this->GetBounds(bounds[0], bounds[1], bounds[2], bounds[3], bounds[4], bounds[5]);
501}
502
504{
505 this->Reset();
506}
507
508inline vtkBoundingBox::vtkBoundingBox(const double bounds[6])
509{
510 this->Reset();
511 this->SetBounds(bounds);
512}
513
515 double xMin, double xMax, double yMin, double yMax, double zMin, double zMax)
516{
517 this->Reset();
518 this->SetBounds(xMin, xMax, yMin, yMax, zMin, zMax);
519}
520
522{
523 this->MinPnt[0] = bbox.MinPnt[0];
524 this->MinPnt[1] = bbox.MinPnt[1];
525 this->MinPnt[2] = bbox.MinPnt[2];
526
527 this->MaxPnt[0] = bbox.MaxPnt[0];
528 this->MaxPnt[1] = bbox.MaxPnt[1];
529 this->MaxPnt[2] = bbox.MaxPnt[2];
530}
531
532inline vtkBoundingBox::vtkBoundingBox(double center[3], double delta)
533{
534 this->Reset();
535 this->AddPoint(center);
536 this->Inflate(delta);
537}
538
540{
541 this->MinPnt[0] = bbox.MinPnt[0];
542 this->MinPnt[1] = bbox.MinPnt[1];
543 this->MinPnt[2] = bbox.MinPnt[2];
544
545 this->MaxPnt[0] = bbox.MaxPnt[0];
546 this->MaxPnt[1] = bbox.MaxPnt[1];
547 this->MaxPnt[2] = bbox.MaxPnt[2];
548 return *this;
549}
550
551inline bool vtkBoundingBox::operator==(const vtkBoundingBox& bbox) const
552{
553 return ((this->MinPnt[0] == bbox.MinPnt[0]) && (this->MinPnt[1] == bbox.MinPnt[1]) &&
554 (this->MinPnt[2] == bbox.MinPnt[2]) && (this->MaxPnt[0] == bbox.MaxPnt[0]) &&
555 (this->MaxPnt[1] == bbox.MaxPnt[1]) && (this->MaxPnt[2] == bbox.MaxPnt[2]));
556}
557
558inline bool vtkBoundingBox::operator!=(const vtkBoundingBox& bbox) const
559{
560 return !((*this) == bbox);
561}
562
563inline void vtkBoundingBox::SetMinPoint(double p[3])
564{
565 this->SetMinPoint(p[0], p[1], p[2]);
566}
567
568inline void vtkBoundingBox::SetMaxPoint(double p[3])
569{
570 this->SetMaxPoint(p[0], p[1], p[2]);
571}
572
573inline void vtkBoundingBox::GetMinPoint(double& x, double& y, double& z) const
574{
575 x = this->MinPnt[0];
576 y = this->MinPnt[1];
577 z = this->MinPnt[2];
578}
579
580inline void vtkBoundingBox::GetMaxPoint(double& x, double& y, double& z) const
581{
582 x = this->MaxPnt[0];
583 y = this->MaxPnt[1];
584 z = this->MaxPnt[2];
585}
586
587inline vtkTypeBool vtkBoundingBox::ContainsPoint(double px, double py, double pz) const
588{
589 if ((px < this->MinPnt[0]) || (px > this->MaxPnt[0]))
590 {
591 return 0;
592 }
593 if ((py < this->MinPnt[1]) || (py > this->MaxPnt[1]))
594 {
595 return 0;
596 }
597 if ((pz < this->MinPnt[2]) || (pz > this->MaxPnt[2]))
598 {
599 return 0;
600 }
601 return 1;
602}
603
604inline vtkTypeBool vtkBoundingBox::ContainsPoint(const double p[3]) const
605{
606 return this->ContainsPoint(p[0], p[1], p[2]);
607}
608
609template <class PointT>
610inline bool vtkBoundingBox::ContainsPoint(const PointT& p) const
611{
612 return this->ContainsPoint(p[0], p[1], p[2]);
613}
614
615inline void vtkBoundingBox::GetCorner(int corner, double p[3]) const
616{
617 if ((corner < 0) || (corner > 7))
618 {
619 p[0] = VTK_DOUBLE_MAX;
620 p[1] = VTK_DOUBLE_MAX;
621 p[2] = VTK_DOUBLE_MAX;
622 return; // out of bounds
623 }
624
625 int ix = (corner & 1) ? 1 : 0; // 0,1,0,1,0,1,0,1
626 int iy = ((corner >> 1) & 1) ? 1 : 0; // 0,0,1,1,0,0,1,1
627 int iz = (corner >> 2) ? 1 : 0; // 0,0,0,0,1,1,1,1
628
629 const double* pts[2] = { this->MinPnt, this->MaxPnt };
630 p[0] = pts[ix][0];
631 p[1] = pts[iy][1];
632 p[2] = pts[iz][2];
633}
634
635VTK_ABI_NAMESPACE_END
636#endif
637// VTK-HeaderTest-Exclude: vtkBoundingBox.h
void ScaleAboutCenter(double s)
Scale each dimension of the box by some given factor, with the origin of the bounding box the center ...
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...
double GetDiagonalLength2() const
Return the length of the diagonal.
static void ComputeBounds(vtkPoints *pts, const long long *ptIds, long long numPointIds, double bounds[6])
Compute the bounding box from an array of vtkPoints.
void Scale(double s[3])
Scale each dimension of the box by some given factor.
void ClampPoint(double point[3])
Clamp point so it is contained inside box.
void AddBounds(const double bounds[6])
Adjust the bounding box so it contains the specified bounds (defined by the VTK representation (xmin,...
vtkIdType ComputeDivisions(vtkIdType totalBins, double bounds[6], int divs[3]) const
Compute the number of divisions in the x-y-z directions given a psoitive, target number of total bins...
int IntersectBox(const vtkBoundingBox &bbox)
Intersect this box with bbox.
const double * GetMinPoint() const
Get the minimum point of the bounding box.
double GetDiagonalLength() const
Return the length of the diagonal.
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.
double GetMaxLength() const
Return the maximum length of the box.
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.
void GetDistance(double point[3], double distance[3])
For each axis, get the minimum distance to put the point inside the box.
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 ...
void Reset()
Returns the box to its initialized state.
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 ComputeBounds(vtkPoints *pts, const long *ptIds, long numPointIds, double bounds[6])
Compute the bounding box from an array of vtkPoints.
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.
static void ClampDivisions(vtkIdType targetBins, int divs[3])
Clamp the number of divisions to be less than or equal to a target number of bins,...
vtkBoundingBox()
Construct a bounding box with the min point set to VTK_DOUBLE_MAX and the max point set to VTK_DOUBLE...
static void ComputeBounds(vtkPoints *pts, const int *ptIds, int numPointIds, double bounds[6])
Compute the bounding box from an array of vtkPoints.
void GetLengths(double lengths[3]) const
Get the length of each side of the box.
void Inflate(double delta)
Expand the bounding box.
void ComputeBounds(vtkPoints *pts, unsigned char *ptUses)
Compute the bounding box from an array of vtkPoints.
void InflateSlice(double delta)
Expand the bounding box.
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 Translate(double motion[3])
Translate box from motion.
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:30
int vtkTypeBool
Definition vtkABI.h:64
void Reset()
Initializes the field list to empty.
virtual double * GetBounds()
Compute the bounding box of the current cell in `bounds' in global coordinates.
bool VTKCOMMONCORE_EXPORT operator==(const std::string &a, const vtkStringToken &b)
bool VTKCOMMONCORE_EXPORT operator!=(const std::string &a, const vtkStringToken &b)
virtual void SetBounds(double, double, double, double, double, double)
int vtkIdType
Definition vtkType.h:332
#define VTK_DOUBLE_MIN
Definition vtkType.h:170
#define VTK_DOUBLE_MAX
Definition vtkType.h:171
int Inflate(double dist) override
Inflates voxel by moving every faces by dist.
#define VTK_SIZEHINT(...)