Loading [MathJax]/extensions/tex2jax.js
Mesh Oriented datABase  (version 5.5.1)
An array-based unstructured mesh library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
moab::OrientedBox Class Reference

Oriented bounding box. More...

#include <OrientedBox.hpp>

+ Collaboration diagram for moab::OrientedBox:

Classes

struct  CovarienceData
 

Public Member Functions

 OrientedBox ()
 
 OrientedBox (const Matrix3 &axes_mat, const CartVect &center)
 
 OrientedBox (const CartVect axes_in[3], const CartVect &center)
 
double inner_radius () const
 radius of inscribed sphere More...
 
double outer_radius () const
 radius of circumscribed sphere More...
 
double outer_radius_squared (const double reps) const
 square of (radius+at least epsilon) of circumsphere More...
 
double inner_radius_squared (const double reps) const
 square of (radius-epsilon) of inscribed sphere More...
 
double volume () const
 volume of box More...
 
CartVect dimensions () const
 number of dimensions for which box is not flat More...
 
double area () const
 largest side area More...
 
CartVect axis (int index) const
 get unit vector in direction of axis More...
 
CartVect scaled_axis (int index) const
 get vector in direction of axis, scaled to its true length More...
 
bool contained (const CartVect &point, double tolerance) const
 
bool intersect_ray (const CartVect &ray_start_point, const CartVect &ray_unit_direction, const double distance_tolerance, const double *nonnegatve_ray_len=0, const double *negative_ray_len=0) const
 
void closest_location_in_box (const CartVect &input_position, CartVect &output_position) const
 Find closest position on/within box to input position. More...
 
ErrorCode make_hex (EntityHandle &hex, Interface *instance)
 Construct a hexahedral element with the same shape as this box. More...
 

Static Public Member Functions

static ErrorCode tag_handle (Tag &handle_out, Interface *instance, const char *name)
 get tag handle for storing oriented box More...
 
static ErrorCode compute_from_vertices (OrientedBox &result, Interface *instance, const Range &vertices)
 Calculate an oriented box from a set of vertices. More...
 
static ErrorCode compute_from_2d_cells (OrientedBox &result, Interface *instance, const Range &elements)
 Calculate an oriented box from a set of 2D elements. More...
 
static ErrorCode covariance_data_from_tris (CovarienceData &result, Interface *moab_instance, const Range &elements)
 
static ErrorCode compute_from_covariance_data (OrientedBox &result, Interface *moab_instance, const CovarienceData *orient_array, unsigned orient_array_length, const Range &vertices)
 
static ErrorCode compute_from_covariance_data (OrientedBox &result, Interface *moab_instance, CovarienceData &orientation_data, const Range &vertices)
 

Public Attributes

CartVect center
 Box center. More...
 
Matrix3 axes
 Box axes, unit vectors sorted by extent of box along axis. More...
 
CartVect length
 distance from center to plane along each axis More...
 
double radius
 outer radius (1/2 diagonal length) of box More...
 

Private Member Functions

void order_axes_by_length (double ax1_len, double ax2_len, double ax3_len)
 orders the box axes by the given lengths for each axis More...
 

Detailed Description

Oriented bounding box.

Definition at line 40 of file OrientedBox.hpp.

Constructor & Destructor Documentation

◆ OrientedBox() [1/3]

moab::OrientedBox::OrientedBox ( )
inline

Definition at line 57 of file OrientedBox.hpp.

57 : radius( 0.0 ) {}

Referenced by compute_from_covariance_data().

◆ OrientedBox() [2/3]

moab::OrientedBox::OrientedBox ( const Matrix3 axes_mat,
const CartVect center 
)

Definition at line 129 of file OrientedBox.cpp.

129  : center( mid ), axes( axes_mat ) 130 { 131  order_axes_by_length( axes.col( 0 ).length(), axes.col( 1 ).length(), axes.col( 2 ).length() ); 132 }

References axes, moab::Matrix3::col(), moab::CartVect::length(), and order_axes_by_length().

◆ OrientedBox() [3/3]

moab::OrientedBox::OrientedBox ( const CartVect  axes_in[3],
const CartVect center 
)

Definition at line 121 of file OrientedBox.cpp.

121  : center( mid ) 122 { 123  124  axes = Matrix3( axes_in[0], axes_in[1], axes_in[2], false ); 125  126  order_axes_by_length( axes_in[0].length(), axes_in[1].length(), axes_in[2].length() ); 127 }

References axes, length, and order_axes_by_length().

Member Function Documentation

◆ area()

double moab::OrientedBox::area ( ) const
inline

largest side area

Definition at line 228 of file OrientedBox.hpp.

229 { 230 #if MB_ORIENTED_BOX_UNIT_VECTORS 231  return 4 * length[1] * length[2]; 232 #else 233  return 4 * ( axes.col( 1 ) * axes.col( 2 ) ).length(); 234 #endif 235 }

References axes, moab::Matrix3::col(), and length.

Referenced by moab::OrientedBoxTreeTool::recursive_stats().

◆ axis()

CartVect moab::OrientedBox::axis ( int  index) const
inline

get unit vector in direction of axis

Definition at line 237 of file OrientedBox.hpp.

238 { 239  return axes.col( index ); 240 }

References axes, and moab::Matrix3::col().

Referenced by moab::TreeNodePrinter::print_geometry().

◆ closest_location_in_box()

void moab::OrientedBox::closest_location_in_box ( const CartVect input_position,
CartVect output_position 
) const

Find closest position on/within box to input position.

Find the closest position in the solid box to the input position. If the input position is on or within the box, then the output position will be the same as the input position. If the input position is outside the box, the outside position will be the closest point on the box boundary to the input position.

Definition at line 658 of file OrientedBox.cpp.

659 { 660  // get coordinates on box axes 661  const CartVect from_center = input_position - center; 662  663 #if MB_ORIENTED_BOX_UNIT_VECTORS 664  CartVect local( from_center % axes.col( 0 ), from_center % axes.col( 1 ), from_center % axes.col( 2 ) ); 665  666  for( int i = 0; i < 3; ++i ) 667  { 668  if( local[i] < -length[i] ) 669  local[i] = -length[i]; 670  else if( local[i] > length[i] ) 671  local[i] = length[i]; 672  } 673 #else 674  CartVect local( ( from_center % axes.col( 0 ) ) / ( axes.col( 0 ) % axes.col( 0 ) ), 675  ( from_center % axes.col( 1 ) ) / ( axes.col( 1 ) % axes.col( 1 ) ), 676  ( from_center % axes.col( 2 ) ) / ( axes.col( 2 ) % axes.col( 2 ) ) ); 677  678  for( int i = 0; i < 3; ++i ) 679  { 680  if( local[i] < -1.0 ) 681  local[i] = -1.0; 682  else if( local[i] > 1.0 ) 683  local[i] = 1.0; 684  } 685 #endif 686  687  output_position = center + local[0] * axes.col( 0 ) + local[1] * axes.col( 1 ) + local[2] * axes.col( 2 ); 688 }

References axes, center, moab::Matrix3::col(), and length.

Referenced by moab::OrientedBoxTreeTool::closest_to_location(), and moab::OrientedBoxTreeTool::sphere_intersect_triangles().

◆ compute_from_2d_cells()

ErrorCode moab::OrientedBox::compute_from_2d_cells ( OrientedBox result,
Interface instance,
const Range elements 
)
static

Calculate an oriented box from a set of 2D elements.

Definition at line 311 of file OrientedBox.cpp.

312 { 313  // Get orientation data from elements 314  CovarienceData data; 315  ErrorCode rval = covariance_data_from_tris( data, instance, elements ); 316  if( MB_SUCCESS != rval ) return rval; 317  318  // get vertices from elements 319  Range points; 320  rval = instance->get_adjacencies( elements, 0, false, points, Interface::UNION ); 321  if( MB_SUCCESS != rval ) return rval; 322  323  // Calculate box given points and orientation data 324  return compute_from_covariance_data( result, instance, data, points ); 325 }

References compute_from_covariance_data(), covariance_data_from_tris(), ErrorCode, moab::Interface::get_adjacencies(), MB_SUCCESS, and moab::Interface::UNION.

Referenced by moab::OrientedBoxTreeTool::build_tree().

◆ compute_from_covariance_data() [1/2]

ErrorCode moab::OrientedBox::compute_from_covariance_data ( OrientedBox result,
Interface moab_instance,
const CovarienceData orient_array,
unsigned  orient_array_length,
const Range vertices 
)
static

Calculate an OrientedBox given an array of CovarienceData and the list of vertices the box is to bound.

Definition at line 372 of file OrientedBox.cpp.

377 { 378  // Sum input CovarienceData structures 379  CovarienceData data_sum( Matrix3( 0.0 ), CartVect( 0.0 ), 0.0 ); 380  for( const CovarienceData* const end = data + data_length; data != end; ++data ) 381  { 382  data_sum.matrix += data->matrix; 383  data_sum.center += data->center; 384  data_sum.area += data->area; 385  } 386  // Compute box from sum of structs 387  return compute_from_covariance_data( result, moab_instance, data_sum, vertices ); 388 }

References moab::OrientedBox::CovarienceData::area, moab::OrientedBox::CovarienceData::center, and moab::OrientedBox::CovarienceData::matrix.

Referenced by moab::OrientedBoxTreeTool::build_sets(), and compute_from_2d_cells().

◆ compute_from_covariance_data() [2/2]

ErrorCode moab::OrientedBox::compute_from_covariance_data ( OrientedBox result,
Interface moab_instance,
CovarienceData orientation_data,
const Range vertices 
)
static

Calculate an OrientedBox given a CovarienceData struct and the list of points the box is to bound.

Definition at line 327 of file OrientedBox.cpp.

331 { 332  if( data.area <= 0.0 ) 333  { 334  Matrix3 empty_axes( 0.0 ); 335  result = OrientedBox( empty_axes, CartVect( 0. ) ); 336  return MB_SUCCESS; 337  } 338  339  // get center from sum 340  result.center = data.center / data.area; 341  342  // get covariance matrix from moments 343  data.matrix /= 12 * data.area; 344  data.matrix -= outer_product( result.center, result.center ); 345  346  // get axes (Eigenvectors) from covariance matrix 347  CartVect lambda; 348  data.matrix.eigen_decomposition( lambda, result.axes ); 349  350  // We now have only the axes. Calculate proper center 351  // and extents for enclosed points. 352  return box_from_axes( result, instance, vertices ); 353 }

References moab::OrientedBox::CovarienceData::area, axes, moab::box_from_axes(), center, moab::OrientedBox::CovarienceData::center, moab::Matrix3::eigen_decomposition(), moab::OrientedBox::CovarienceData::matrix, MB_SUCCESS, OrientedBox(), and moab::outer_product().

◆ compute_from_vertices()

ErrorCode moab::OrientedBox::compute_from_vertices ( OrientedBox result,
Interface instance,
const Range vertices 
)
static

Calculate an oriented box from a set of vertices.

Definition at line 231 of file OrientedBox.cpp.

232 { 233  const Range::iterator begin = vertices.lower_bound( MBVERTEX ); 234  const Range::iterator end = vertices.upper_bound( MBVERTEX ); 235  size_t count = 0; 236  237  // compute mean 238  CartVect v; 239  result.center = CartVect( 0, 0, 0 ); 240  for( Range::iterator i = begin; i != end; ++i ) 241  { 242  ErrorCode rval = instance->get_coords( &*i, 1, v.array() ); 243  if( MB_SUCCESS != rval ) return rval; 244  result.center += v; 245  ++count; 246  } 247  result.center /= count; 248  249  // compute covariance matrix 250  Matrix3 a( 0.0 ); 251  for( Range::iterator i = begin; i != end; ++i ) 252  { 253  ErrorCode rval = instance->get_coords( &*i, 1, v.array() ); 254  if( MB_SUCCESS != rval ) return rval; 255  256  v -= result.center; 257  a += outer_product( v, v ); 258  } 259  a /= count; 260  261  // Get axes (Eigenvectors) from covariance matrix 262  CartVect lambda; 263  a.eigen_decomposition( lambda, result.axes ); 264  265  // Calculate center and extents of box given orientation defined by axes 266  return box_from_axes( result, instance, vertices ); 267 }

References moab::CartVect::array(), axes, moab::box_from_axes(), center, moab::Matrix3::eigen_decomposition(), ErrorCode, moab::Interface::get_coords(), moab::Range::lower_bound(), MB_SUCCESS, MBVERTEX, moab::outer_product(), and moab::Range::upper_bound().

◆ contained()

bool moab::OrientedBox::contained ( const CartVect point,
double  tolerance 
) const

Test if point is contained in box

Definition at line 355 of file OrientedBox.cpp.

356 { 357  CartVect from_center = point - center; 358 #if MB_ORIENTED_BOX_UNIT_VECTORS 359  return fabs( from_center % axes.col( 0 ) ) - length[0] <= tol && 360  fabs( from_center % axes.col( 1 ) ) - length[1] <= tol && 361  fabs( from_center % axes.col( 2 ) ) - length[2] <= tol; 362 #else 363  for( int i = 0; i < 3; ++i ) 364  { 365  double length = axes.col( i ).length(); 366  if( fabs( from_center % axes.col( i ) ) - length * length > length * tol ) return false; 367  } 368  return true; 369 #endif 370 }

References axes, center, moab::Matrix3::col(), moab::CartVect::length(), and length.

Referenced by TriCounter::visit().

◆ covariance_data_from_tris()

ErrorCode moab::OrientedBox::covariance_data_from_tris ( CovarienceData result,
Interface moab_instance,
const Range elements 
)
static

Calculate a CovarienceData struct from a list of triangles

Definition at line 269 of file OrientedBox.cpp.

270 { 271  ErrorCode rval; 272  const Range::iterator begin = elements.lower_bound( CN::TypeDimensionMap[2].first ); 273  const Range::iterator end = elements.lower_bound( CN::TypeDimensionMap[3].first ); 274  275  // compute mean and moments 276  result.matrix = Matrix3( 0.0 ); 277  result.center = CartVect( 0.0 ); 278  result.area = 0.0; 279  for( Range::iterator i = begin; i != end; ++i ) 280  { 281  const EntityHandle* conn = NULL; 282  int conn_len = 0; 283  rval = instance->get_connectivity( *i, conn, conn_len ); 284  if( MB_SUCCESS != rval ) return rval; 285  286  // for each triangle in the 2-D cell 287  for( int j = 2; j < conn_len; ++j ) 288  { 289  EntityHandle vertices[3] = { conn[0], conn[j - 1], conn[j] }; 290  CartVect coords[3]; 291  rval = instance->get_coords( vertices, 3, coords[0].array() ); 292  if( MB_SUCCESS != rval ) return rval; 293  294  // edge vectors 295  const CartVect edge0 = coords[1] - coords[0]; 296  const CartVect edge1 = coords[2] - coords[0]; 297  const CartVect centroid = ( coords[0] + coords[1] + coords[2] ) / 3; 298  const double tri_area2 = ( edge0 * edge1 ).length(); 299  result.area += tri_area2; 300  result.center += tri_area2 * centroid; 301  302  result.matrix += 303  tri_area2 * ( 9 * outer_product( centroid, centroid ) + outer_product( coords[0], coords[0] ) + 304  outer_product( coords[1], coords[1] ) + outer_product( coords[2], coords[2] ) ); 305  } // for each triangle 306  } // for each element 307  308  return MB_SUCCESS; 309 }

References moab::OrientedBox::CovarienceData::area, moab::OrientedBox::CovarienceData::center, ErrorCode, moab::GeomUtil::first(), moab::Interface::get_connectivity(), moab::Interface::get_coords(), length, moab::Range::lower_bound(), moab::OrientedBox::CovarienceData::matrix, MB_SUCCESS, moab::outer_product(), and moab::CN::TypeDimensionMap.

Referenced by compute_from_2d_cells(), and moab::OrientedBoxTreeTool::join_trees().

◆ dimensions()

CartVect moab::OrientedBox::dimensions ( ) const
inline

number of dimensions for which box is not flat

Definition at line 219 of file OrientedBox.hpp.

220 { 221 #if MB_ORIENTED_BOX_UNIT_VECTORS 222  return 2.0 * length; 223 #else 224  return 2.0 * CartVect( axes.col( 0 ).length(), axes.col( 1 ).length(), axes.col( 2 ).length() ); 225 #endif 226 }

References axes, moab::Matrix3::col(), moab::CartVect::length(), and length.

Referenced by moab::TreeNodePrinter::print_geometry(), and moab::OrientedBoxTreeTool::recursive_stats().

◆ inner_radius()

double moab::OrientedBox::inner_radius ( ) const
inline

radius of inscribed sphere

Definition at line 163 of file OrientedBox.hpp.

164 { 165 #if MB_ORIENTED_BOX_UNIT_VECTORS 166  return length[0]; 167 #else 168  return axes.col( 0 ).length(); 169 #endif 170 }

References axes, moab::Matrix3::col(), moab::CartVect::length(), and length.

Referenced by moab::TreeNodePrinter::print_geometry(), and moab::OrientedBoxTreeTool::recursive_stats().

◆ inner_radius_squared()

double moab::OrientedBox::inner_radius_squared ( const double  reps) const
inline

square of (radius-epsilon) of inscribed sphere

Definition at line 199 of file OrientedBox.hpp.

200 { 201 #if MB_ORIENTED_BOX_UNIT_VECTORS 202  return ( length[0] - reps ) * ( length[0] - reps ); 203 #else 204  CartVect tmp = axes.col( 0 ); 205  tmp -= CartVect( reps, reps, reps ); 206  return ( tmp % tmp ); 207 #endif 208 }

References axes, moab::Matrix3::col(), and length.

Referenced by intersect_ray().

◆ intersect_ray()

bool moab::OrientedBox::intersect_ray ( const CartVect ray_start_point,
const CartVect ray_unit_direction,
const double  distance_tolerance,
const double *  nonnegatve_ray_len = 0,
const double *  negative_ray_len = 0 
) const

Test for intersection of a ray (or line segment) with this box. Ray length limits are used to optimize Monte Carlo particle tracking.

Parameters
ray_start_pointThe base point of the ray
ray_unit_directionThe direction of the ray (must be unit length)
distance_toleranceTolerance to use in intersection checks
nonnegative_ray_lenOptional length of ray in forward direction
negative_ray_lenOptional length of ray in reverse direction

Definition at line 498 of file OrientedBox.cpp.

503 { 504  // test distance from box center to line 505  const CartVect cx = center - ray_origin; 506  const double dist_s = cx % ray_direction; 507  const double dist_sq = cx % cx - ( dist_s * dist_s ); 508  const double max_diagsq = outer_radius_squared( reps ); 509  510  // For the largest sphere, no intersections exist if discriminant is negative. 511  // Geometrically, if distance from box center to line is greater than the 512  // longest diagonal, there is no intersection. 513  // manipulate the discriminant: 0 > dist_s*dist_s - cx%cx + max_diagsq 514  if( dist_sq > max_diagsq ) return false; 515  516  // If the closest possible intersection must be closer than nonneg_ray_len. Be 517  // careful with absolute value, squaring distances, and subtracting squared 518  // distances. 519  if( nonneg_ray_len ) 520  { 521  assert( 0 <= *nonneg_ray_len ); 522  double max_len; 523  if( neg_ray_len ) 524  { 525  assert( 0 >= *neg_ray_len ); 526  max_len = std::max( *nonneg_ray_len, -( *neg_ray_len ) ); 527  } 528  else 529  { 530  max_len = *nonneg_ray_len; 531  } 532  const double temp = fabs( dist_s ) - max_len; 533  if( 0.0 < temp && temp * temp > max_diagsq ) return false; 534  } 535  536  // if smaller than shortest diagonal, we do hit 537  if( dist_sq < inner_radius_squared( reps ) ) 538  { 539  // nonnegative direction 540  if( dist_s >= 0.0 ) 541  { 542  if( nonneg_ray_len ) 543  { 544  if( *nonneg_ray_len > dist_s ) return true; 545  } 546  else 547  { 548  return true; 549  } 550  // negative direction 551  } 552  else 553  { 554  if( neg_ray_len && *neg_ray_len < dist_s ) return true; 555  } 556  } 557  558  // get transpose of axes 559  Matrix3 B = axes.transpose(); 560  561  // transform ray to box coordintae system 562  CartVect par_pos = B * ( ray_origin - center ); 563  CartVect par_dir = B * ray_direction; 564  565  // Fast Rejection Test: Ray will not intersect if it is going away from the box. 566  // This will not work for rays with neg_ray_len. length[0] is half of box width 567  // along axes.col(0). 568  const double half_x = length[0] + reps; 569  const double half_y = length[1] + reps; 570  const double half_z = length[2] + reps; 571  if( !neg_ray_len ) 572  { 573  if( ( par_pos[0] > half_x && par_dir[0] >= 0 ) || ( par_pos[0] < -half_x && par_dir[0] <= 0 ) ) return false; 574  575  if( ( par_pos[1] > half_y && par_dir[1] >= 0 ) || ( par_pos[1] < -half_y && par_dir[1] <= 0 ) ) return false; 576  577  if( ( par_pos[2] > half_z && par_dir[2] >= 0 ) || ( par_pos[2] < -half_z && par_dir[2] <= 0 ) ) return false; 578  } 579  580  // test if ray_origin is inside box 581  if( par_pos[0] <= half_x && par_pos[0] >= -half_x && par_pos[1] <= half_y && par_pos[1] >= -half_y && 582  par_pos[2] <= half_z && par_pos[2] >= -half_z ) 583  return true; 584  585  // test two xy plane 586  if( fabs( par_dir[0] * ( half_z - par_pos[2] ) + par_dir[2] * par_pos[0] ) <= 587  fabs( par_dir[2] * half_x ) && // test against x extents using z 588  fabs( par_dir[1] * ( half_z - par_pos[2] ) + par_dir[2] * par_pos[1] ) <= 589  fabs( par_dir[2] * half_y ) && // test against y extents using z 590  check_ray_limits( par_pos[2], par_dir[2], half_z, nonneg_ray_len, neg_ray_len ) ) 591  return true; 592  if( fabs( par_dir[0] * ( -half_z - par_pos[2] ) + par_dir[2] * par_pos[0] ) <= fabs( par_dir[2] * half_x ) && 593  fabs( par_dir[1] * ( -half_z - par_pos[2] ) + par_dir[2] * par_pos[1] ) <= fabs( par_dir[2] * half_y ) && 594  check_ray_limits( par_pos[2], par_dir[2], -half_z, nonneg_ray_len, neg_ray_len ) ) 595  return true; 596  597  // test two xz plane 598  if( fabs( par_dir[0] * ( half_y - par_pos[1] ) + par_dir[1] * par_pos[0] ) <= fabs( par_dir[1] * half_x ) && 599  fabs( par_dir[2] * ( half_y - par_pos[1] ) + par_dir[1] * par_pos[2] ) <= fabs( par_dir[1] * half_z ) && 600  check_ray_limits( par_pos[1], par_dir[1], half_y, nonneg_ray_len, neg_ray_len ) ) 601  return true; 602  if( fabs( par_dir[0] * ( -half_y - par_pos[1] ) + par_dir[1] * par_pos[0] ) <= fabs( par_dir[1] * half_x ) && 603  fabs( par_dir[2] * ( -half_y - par_pos[1] ) + par_dir[1] * par_pos[2] ) <= fabs( par_dir[1] * half_z ) && 604  check_ray_limits( par_pos[1], par_dir[1], -half_y, nonneg_ray_len, neg_ray_len ) ) 605  return true; 606  607  // test two yz plane 608  if( fabs( par_dir[1] * ( half_x - par_pos[0] ) + par_dir[0] * par_pos[1] ) <= fabs( par_dir[0] * half_y ) && 609  fabs( par_dir[2] * ( half_x - par_pos[0] ) + par_dir[0] * par_pos[2] ) <= fabs( par_dir[0] * half_z ) && 610  check_ray_limits( par_pos[0], par_dir[0], half_x, nonneg_ray_len, neg_ray_len ) ) 611  return true; 612  if( fabs( par_dir[1] * ( -half_x - par_pos[0] ) + par_dir[0] * par_pos[1] ) <= fabs( par_dir[0] * half_y ) && 613  fabs( par_dir[2] * ( -half_x - par_pos[0] ) + par_dir[0] * par_pos[2] ) <= fabs( par_dir[0] * half_z ) && 614  check_ray_limits( par_pos[0], par_dir[0], -half_x, nonneg_ray_len, neg_ray_len ) ) 615  return true; 616  617  return false; 618 }

References axes, center, moab::check_ray_limits(), inner_radius_squared(), length, outer_radius_squared(), and moab::Matrix3::transpose().

◆ make_hex()

ErrorCode moab::OrientedBox::make_hex ( EntityHandle hex,
Interface instance 
)

Construct a hexahedral element with the same shape as this box.

Definition at line 620 of file OrientedBox.cpp.

621 { 622  ErrorCode rval; 623  int signs[8][3] = { { -1, -1, -1 }, { 1, -1, -1 }, { 1, 1, -1 }, { -1, 1, -1 }, 624  { -1, -1, 1 }, { 1, -1, 1 }, { 1, 1, 1 }, { -1, 1, 1 } }; 625  626  std::vector< EntityHandle > vertices; 627  for( int i = 0; i < 8; ++i ) 628  { 629  CartVect coords( center ); 630  for( int j = 0; j < 3; ++j ) 631  { 632 #if MB_ORIENTED_BOX_UNIT_VECTORS 633  coords += signs[i][j] * ( axes.col( j ) * length[j] ); 634 #else 635  coords += signs[i][j] * axes.col( j ); 636 #endif 637  } 638  EntityHandle handle; 639  rval = instance->create_vertex( coords.array(), handle ); 640  if( MB_SUCCESS != rval ) 641  { 642  instance->delete_entities( &vertices[0], vertices.size() ); 643  return rval; 644  } 645  vertices.push_back( handle ); 646  } 647  648  rval = instance->create_element( MBHEX, &vertices[0], vertices.size(), hex ); 649  if( MB_SUCCESS != rval ) 650  { 651  instance->delete_entities( &vertices[0], vertices.size() ); 652  return rval; 653  } 654  655  return MB_SUCCESS; 656 }

References moab::CartVect::array(), axes, center, moab::Matrix3::col(), moab::Interface::create_element(), moab::Interface::create_vertex(), moab::Interface::delete_entities(), ErrorCode, length, MB_SUCCESS, and MBHEX.

◆ order_axes_by_length()

void moab::OrientedBox::order_axes_by_length ( double  ax1_len,
double  ax2_len,
double  ax3_len 
)
private

orders the box axes by the given lengths for each axis

Definition at line 85 of file OrientedBox.cpp.

86 { 87  88  CartVect len( ax1_len, ax2_len, ax3_len ); 89  90  if( len[2] < len[1] ) 91  { 92  if( len[2] < len[0] ) 93  { 94  std::swap( len[0], len[2] ); 95  axes.swapcol( 0, 2 ); 96  } 97  } 98  else if( len[1] < len[0] ) 99  { 100  std::swap( len[0], len[1] ); 101  axes.swapcol( 0, 1 ); 102  } 103  if( len[1] > len[2] ) 104  { 105  std::swap( len[1], len[2] ); 106  axes.swapcol( 1, 2 ); 107  } 108  109 #if MB_ORIENTED_BOX_UNIT_VECTORS 110  length = len; 111  if( len[0] > 0.0 ) axes.colscale( 0, 1.0 / len[0] ); 112  if( len[1] > 0.0 ) axes.colscale( 1, 1.0 / len[1] ); 113  if( len[2] > 0.0 ) axes.colscale( 2, 1.0 / len[2] ); 114 #endif 115  116 #if MB_ORIENTED_BOX_OUTER_RADIUS 117  radius = len.length(); 118 #endif 119 }

References axes, moab::Matrix3::colscale(), moab::CartVect::length(), length, radius, and moab::Matrix3::swapcol().

Referenced by OrientedBox().

◆ outer_radius()

double moab::OrientedBox::outer_radius ( ) const
inline

radius of circumscribed sphere

Definition at line 172 of file OrientedBox.hpp.

173 { 174 #if MB_ORIENTED_BOX_OUTER_RADIUS 175  return radius; 176 #elif MB_ORIENTED_BOX_UNIT_VECTORS 177  return length.length(); 178 #else 179  return ( axes.col( 0 ) + axes.col( 1 ) + axes.col( 2 ) ).length(); 180 #endif 181 }

References axes, moab::Matrix3::col(), moab::CartVect::length(), length, and radius.

Referenced by moab::TreeNodePrinter::print_geometry(), and moab::OrientedBoxTreeTool::recursive_stats().

◆ outer_radius_squared()

double moab::OrientedBox::outer_radius_squared ( const double  reps) const
inline

square of (radius+at least epsilon) of circumsphere

Definition at line 184 of file OrientedBox.hpp.

185 { 186 #if MB_ORIENTED_BOX_OUTER_RADIUS 187  return ( radius + reps ) * ( radius + reps ); 188 #elif MB_ORIENTED_BOX_UNIT_VECTORS 189  CartVect tmp( length[0] + reps, length[1] + reps, length[2] + reps ); 190  return tmp % tmp; 191 #else 192  CartVect half_diag = axes.col( 0 ) + axes.col( 1 ) + axes.col( 2 ); 193  half_diag += CartVect( reps, reps, reps ); 194  return half_diag % half_diag; 195 #endif 196 }

References axes, moab::Matrix3::col(), length, and radius.

Referenced by intersect_ray().

◆ scaled_axis()

CartVect moab::OrientedBox::scaled_axis ( int  index) const
inline

get vector in direction of axis, scaled to its true length

Definition at line 242 of file OrientedBox.hpp.

243 { 244 #if MB_ORIENTED_BOX_UNIT_VECTORS 245  return length[index] * axes.col( index ); 246 #else 247  return axes.col( index ); 248 #endif 249 }

References axes, moab::Matrix3::col(), and length.

Referenced by moab::OrientedBoxTreeTool::box().

◆ tag_handle()

ErrorCode moab::OrientedBox::tag_handle ( Tag handle_out,
Interface instance,
const char *  name 
)
static

get tag handle for storing oriented box

Get the handle for the tag with the specified name and check that the tag is appropriate for storing instances of OrientedBox. The resulting tag may be used to store instances of OrientedBox directly.

Parameters
handle_outThe TagHandle, passed back to caller
nameThe tag name
createIf true, tag will be created if it does not exist

Definition at line 134 of file OrientedBox.cpp.

135 { 136  // We're going to assume this when mapping the OrientedBox 137  // to tag data, so assert it. 138 #if MB_ORIENTED_BOX_OUTER_RADIUS 139  const int rad_size = 1; 140 #else 141  const int rad_size = 0; 142 #endif 143 #if MB_ORIENTED_BOX_UNIT_VECTORS 144  const int SIZE = rad_size + 15; 145 #else 146  const int SIZE = rad_size + 12; 147 #endif 148  assert( sizeof( OrientedBox ) == SIZE * sizeof( double ) ); 149  150  return instance->tag_get_handle( name, SIZE, MB_TYPE_DOUBLE, handle_out, MB_TAG_DENSE | MB_TAG_CREAT ); 151 }

References MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_DOUBLE, and moab::Interface::tag_get_handle().

Referenced by moab::OrientedBoxTreeTool::OrientedBoxTreeTool().

◆ volume()

double moab::OrientedBox::volume ( ) const
inline

volume of box

Definition at line 210 of file OrientedBox.hpp.

211 { 212 #if MB_ORIENTED_BOX_UNIT_VECTORS 213  return 8 * length[0] * length[1] * length[2]; 214 #else 215  return fabs( 8 * axes.col( 0 ) % ( axes.col( 1 ) * axes.col( 2 ) ) ); 216 #endif 217 }

References axes, moab::Matrix3::col(), and length.

Referenced by TriStats::leaf(), moab::OrientedBoxTreeTool::recursive_stats(), and TriStats::TriStats().

Member Data Documentation

◆ axes

◆ center

◆ length

◆ radius

double moab::OrientedBox::radius

outer radius (1/2 diagonal length) of box

Definition at line 54 of file OrientedBox.hpp.

Referenced by moab::box_from_axes(), order_axes_by_length(), outer_radius(), and outer_radius_squared().


The documentation for this class was generated from the following files: