cgma
PointGridSearch Class Reference

#include <PointGridSearch.hpp>

List of all members.

Public Member Functions

 PointGridSearch (DLIList< CubitPoint * > *point_list, DLIList< CubitFacet * > *facet_list, float grid_scale=2.0)
 PointGridSearch (DLIList< CubitPoint * > &point_list, double grid_cell_size, double grid_scale)
virtual ~PointGridSearch ()
CubitVector grid_range_minimum ()
CubitVector grid_range_maximum ()
CubitVector grid_cell_width ()
int number_grid_cells_x ()
int number_grid_cells_y ()
int number_grid_cells_z ()
int number_grid_cells ()
double fraction_empty_cells ()
double average_points_per_occupied_cell ()
int total_points_in_grid_cells ()
CubitVector bounding_range_minimum ()
CubitVector bounding_range_maximum ()
int bounding_cell_minimum_x ()
int bounding_cell_minimum_y ()
int bounding_cell_minimum_z ()
int bounding_cell_maximum_x ()
int bounding_cell_maximum_y ()
int bounding_cell_maximum_z ()
void add_point (CubitPoint *point)
void remove_point (CubitPoint *point)
int in_grid (const CubitVector &position)
void change_cell (CubitPoint *point)
void set_neighborhood_bounds (CubitVector &vec)
void set_neighborhood_bounds (const CubitVector &center, double size)
void set_neighborhood_bounds (CubitVector &center, double x, double y, double z)
void set_neighborhood_bounds (CubitVector &vec_1, CubitVector &vec_2)
void set_neighborhood_bounds (CubitPoint *point, CubitVector &vec)
void set_neighborhood_bounds (CubitPoint *point)
void set_neighborhood_bounds (CubitPoint *point_1, CubitPoint *point_2)
void set_neighborhood_bounds (DLIList< CubitPoint * > &point_list)
void set_neighborhood_bounds (CubitFacet *facet)
void set_neighborhood_bounds (DLIList< CubitFacet * > &facet_list)
void set_neighborhood_bounds_max_edge (const CubitVector &center, double factor=1.0)
void get_neighborhood_points (DLIList< CubitPoint * > &point_list)
void get_neighborhood_facets (DLIList< CubitFacet * > &facet_list)
void get_neighborhood_points_sorted (DLIList< CubitPoint * > &point_list, const CubitVector &center, double cut_off)

Protected Member Functions

int grid_cell_index (CubitPoint *point)
virtual void cell_from_range ()
void bounding_range (DLIList< CubitPoint * > &point_list)
void bounding_range (DLIList< CubitFacet * > &facet_list)
void bounding_range (CubitPoint *point)
void bounding_range (const CubitVector &vec)

Protected Attributes

CubitVector gridRangeMinimum
CubitVector gridRangeMaximum
CubitVector gridCellWidth
CubitVector gridCellWidthInverse
int numberGridCellsX
int numberGridCellsY
int numberGridCellsZ
int numberGridCells
CubitVector boundingRangeMinimum
CubitVector boundingRangeMaximum
int boundingCellMinimumX
int boundingCellMinimumY
int boundingCellMinimumZ
int boundingCellMaximumX
int boundingCellMaximumY
int boundingCellMaximumZ

Private Member Functions

void add_point_to_cell (CubitPoint *point, int i)
void remove_point_from_cell (CubitPoint *point, int i)

Private Attributes

DLIList< CubitPoint * > ** neighborhoodList
double maxEdgeLength

Detailed Description

Definition at line 45 of file PointGridSearch.hpp.


Constructor & Destructor Documentation

PointGridSearch::PointGridSearch ( DLIList< CubitPoint * > *  point_list,
DLIList< CubitFacet * > *  facet_list,
float  grid_scale = 2.0 
)

Definition at line 33 of file PointGridSearch.cpp.

{
    // find the geometric bounding range mesh
  assert( point_list && point_list->size() );
  bounding_range( *point_list );

    //find an estimate of the cell size based on the average square distance
    //of the edges on the surface of the geometry

  double grid_cell_size = 1.;
  
    // edge lengths
  double sum = 0.0;
  int i;
  maxEdgeLength = CUBIT_DBL_MIN;
  assert( facet_list && facet_list->size());
  for ( i = 0; i < facet_list->size(); i++)
  {
    CubitPoint *p1, *p2;
    facet_list->next(i)->get_edge_1(p1, p2);
    CubitVector temp = p2->coordinates() - p1->coordinates();
    double tmp_lsqrd =  temp.length_squared();
    sum += tmp_lsqrd;
  }
  grid_cell_size = sqrt( sum / facet_list->size() );

  maxEdgeLength = grid_cell_size;
    // evaluate grid parameters
  CubitVector cell(1.0, 1.0, 1.0);
  double grid_cell_width = grid_scale * grid_cell_size;
  cell *= ( GRID_EXTENSION + 5.0 ) * grid_cell_width;
  gridRangeMinimum = boundingRangeMinimum - cell;
  gridRangeMaximum = boundingRangeMaximum + cell;
  
  CubitVector range_width = gridRangeMaximum - gridRangeMinimum;
  
  numberGridCellsX = (int) ceil(range_width.x() / grid_cell_width + 0.5);
  numberGridCellsY = (int) ceil(range_width.y() / grid_cell_width + 0.5);
  numberGridCellsZ = (int) ceil(range_width.z() / grid_cell_width + 0.5);
  numberGridCells  = numberGridCellsX * numberGridCellsY * numberGridCellsZ;
  
  gridCellWidth.x(range_width.x() / ((double) numberGridCellsX));
  gridCellWidth.y(range_width.y() / ((double) numberGridCellsY));
  gridCellWidth.z(range_width.z() / ((double) numberGridCellsZ));
  
  gridCellWidthInverse.x(1.0 / gridCellWidth.x());
  gridCellWidthInverse.y(1.0 / gridCellWidth.y());
  gridCellWidthInverse.z(1.0 / gridCellWidth.z());
  
    // allocate neighborhood list array
  neighborhoodList = new DLIList<CubitPoint*>* [numberGridCells];
  assert(neighborhoodList != NULL);
  for ( i = 0; i < numberGridCells; i++)
  {
    neighborhoodList[i] = NULL;
  }
}
PointGridSearch::PointGridSearch ( DLIList< CubitPoint * > &  point_list,
double  grid_cell_size,
double  grid_scale 
)

Definition at line 92 of file PointGridSearch.cpp.

{
    // find the bounding range of the reference geometry
  bounding_range(point_list);
  if ( grid_cell_size < CUBIT_RESABS )
    grid_cell_size = 
      (boundingRangeMaximum - boundingRangeMinimum).length() / 5.;    
  
  CubitVector cell(1.0, 1.0, 1.0);
  double grid_cell_width = grid_scale * grid_cell_size;
  cell *= grid_cell_width;
  gridRangeMinimum = boundingRangeMinimum - cell;
  gridRangeMaximum = boundingRangeMaximum + cell;

  CubitVector range_width = gridRangeMaximum - gridRangeMinimum;

  numberGridCellsX = (int) ceil(range_width.x() / grid_cell_width + 0.5);
  numberGridCellsY = (int) ceil(range_width.y() / grid_cell_width + 0.5);
  numberGridCellsZ = (int) ceil(range_width.z() / grid_cell_width + 0.5);
  numberGridCells  = numberGridCellsX * numberGridCellsY * numberGridCellsZ;

  gridCellWidth.x(range_width.x() / ((double) numberGridCellsX));
  gridCellWidth.y(range_width.y() / ((double) numberGridCellsY));
  gridCellWidth.z(range_width.z() / ((double) numberGridCellsZ));

  gridCellWidthInverse.x(1.0 / gridCellWidth.x());
  gridCellWidthInverse.y(1.0 / gridCellWidth.y());
  gridCellWidthInverse.z(1.0 / gridCellWidth.z());

    // allocate neighborhood list array

  neighborhoodList = new DLIList<CubitPoint*>* [numberGridCells];
  assert(neighborhoodList != NULL);
  int i;
  for ( i = 0; i < numberGridCells; i++)
  {
    neighborhoodList[i] = NULL;
  }
}

Definition at line 133 of file PointGridSearch.cpp.

{
    // delete each DLIList<CubitPoint*> and the array of pointers containing them

    // first, delete the tdcellindex for each point
  int i;
  for ( i = 0; i < numberGridCells; i++) {
    if (!neighborhoodList[i]) continue;
    for (int j = neighborhoodList[i]->size(); j > 0; j--)
      neighborhoodList[i]->get_and_step()->delete_TD(&TDCellIndex::is_cell_index);
  }

    for (i = 0; i < numberGridCells; i++)
      delete neighborhoodList[i];

    delete [] neighborhoodList;
}

Member Function Documentation

Definition at line 290 of file PointGridSearch.cpp.

{
    // append point to the appropriate list

    int i = grid_cell_index(point);
    if (i < 0 || i > numberGridCells){
      return;
    }

    ToolData* td_temp = point->get_TD(&TDCellIndex::is_cell_index);
    TDCellIndex* td_index;
    
    if(td_temp == NULL){
      td_index = new TDCellIndex;
      assert(td_index != NULL);
      point->add_TD(td_index);
      td_index->cell_index(i); 
      add_point_to_cell(point, i);
    }
    else{  // the tool data already exists, avoid multiple listings of point
      td_index = CAST_TO(td_temp, TDCellIndex);
      assert( td_index != NULL);
      int iold = td_index->cell_index();
      
      if(iold != i){
    remove_point_from_cell(point,iold);
    td_index->cell_index(i);
    add_point_to_cell(point, i);
      }
      else{  // iold==i, so check whether the point is already on the list
    if (neighborhoodList[i] && neighborhoodList[i]->move_to(point)) return;
    else add_point_to_cell(point, i);
      }
    }
}
void PointGridSearch::add_point_to_cell ( CubitPoint point,
int  i 
) [private]

Definition at line 667 of file PointGridSearch.cpp.

{
  if (neighborhoodList[i])
  {
    neighborhoodList[i]->append(point);
  }
  else
  {
    neighborhoodList[i] = new DLIList<CubitPoint*>;
    neighborhoodList[i]->append(point);
  }
}

Definition at line 634 of file PointGridSearch.cpp.

{
  int total_points    = 0;
  int occupied_cells = 0;

  for (int i = 0; i < numberGridCells; i++)
  {
    if (neighborhoodList[i])
    {
      total_points += neighborhoodList[i]->size();
      occupied_cells++;
    }
  }

  return ((double) (total_points) / (double) (occupied_cells));
}

Definition at line 205 of file PointGridSearch.hpp.

Definition at line 207 of file PointGridSearch.hpp.

Definition at line 209 of file PointGridSearch.hpp.

Definition at line 199 of file PointGridSearch.hpp.

Definition at line 201 of file PointGridSearch.hpp.

Definition at line 203 of file PointGridSearch.hpp.

void PointGridSearch::bounding_range ( DLIList< CubitPoint * > &  point_list) [protected]

Definition at line 182 of file PointGridSearch.cpp.

{

  if ( !point_list.size() )
    return;
  
    // initialize min and max range values to the first point coordinates
  boundingRangeMinimum = point_list.get()->coordinates();
  boundingRangeMaximum = point_list.get_and_step()->coordinates();
  
    // find the min and max coordinates that completely encloses the
    // point list
  
  for (int i = 1; i < point_list.size(); i++)
    bounding_range(point_list.get_and_step());

}
void PointGridSearch::bounding_range ( DLIList< CubitFacet * > &  facet_list) [protected]

Definition at line 199 of file PointGridSearch.cpp.

{
  int i, j;
  DLIList<CubitPoint*> point_list;
  for ( i = 0; i < facet_list.size(); i++)
  {
    CubitFacet* facet = facet_list.get_and_step();
    
    DLIList<CubitPoint*> temp_point_list;
    facet->points(temp_point_list);
    for (j = 0; j < temp_point_list.size(); j++)
    {
      CubitPoint* point = temp_point_list.get_and_step();
      
      if (!point->marked())
      {
        point->marked(CUBIT_TRUE);
        point_list.append(point);
      }
    }
  }
    // unmark the found points
  
  for (i = 0; i < point_list.size(); i++)
  {
    point_list.get_and_step()->marked(CUBIT_FALSE);
  }
  if ( !point_list.size() )
    return;
  
    // initialize min and max range values to the first point coordinates
  boundingRangeMinimum = point_list.get()->coordinates();
  boundingRangeMaximum = point_list.get_and_step()->coordinates();
  
    // find the min and max coordinates that completely encloses the
    // point list
  
  for (i = 1; i < point_list.size(); i++)
    bounding_range(point_list.get_and_step());

}
void PointGridSearch::bounding_range ( CubitPoint point) [inline, protected]

Definition at line 242 of file PointGridSearch.hpp.

{
  double xpos = point->x();
  double ypos = point->y();
  double zpos = point->z();
  
  if (xpos < boundingRangeMinimum.x())
    boundingRangeMinimum.x(xpos);
  else if (xpos > boundingRangeMaximum.x())
    boundingRangeMaximum.x(xpos);

  if (ypos < boundingRangeMinimum.y())
    boundingRangeMinimum.y(ypos);
  else if (ypos > boundingRangeMaximum.y())
    boundingRangeMaximum.y(ypos);

  if (zpos < boundingRangeMinimum.z())
    boundingRangeMinimum.z(zpos);
  else if (zpos > boundingRangeMaximum.z())
    boundingRangeMaximum.z(zpos);
}
void PointGridSearch::bounding_range ( const CubitVector vec) [inline, protected]

Definition at line 224 of file PointGridSearch.hpp.

{
  if (vec.x() < boundingRangeMinimum.x())
    boundingRangeMinimum.x(vec.x());
  else if (vec.x() > boundingRangeMaximum.x())
    boundingRangeMaximum.x(vec.x());

  if (vec.y() < boundingRangeMinimum.y())
    boundingRangeMinimum.y(vec.y());
  else if (vec.y() > boundingRangeMaximum.y())
    boundingRangeMaximum.y(vec.y());

  if (vec.z() < boundingRangeMinimum.z())
    boundingRangeMinimum.z(vec.z());
  else if (vec.z() > boundingRangeMaximum.z())
    boundingRangeMaximum.z(vec.z());
}

Definition at line 245 of file PointGridSearch.cpp.

{
  int old_cell, current_cell;
  // if cell index has changed then swap lists
  ToolData* td_temp = point->get_TD( &TDCellIndex::is_cell_index);
  TDCellIndex* td_index = CAST_TO(td_temp, TDCellIndex);
  if (td_temp == NULL) {
    td_index = new TDCellIndex;
    assert( td_index != NULL );
    point->add_TD(td_index);
  }
  old_cell = td_index->cell_index();
  
  current_cell = grid_cell_index(point);
  if (current_cell >= 0 && current_cell < numberGridCells)
  {
    if (current_cell != old_cell)
    {
      remove_point_from_cell(point, old_cell);
      add_point_to_cell(point, current_cell);
      td_index->cell_index(current_cell);
    }
  }
}

Definition at line 622 of file PointGridSearch.cpp.

{
  int empty_cell = 0;

  for (int i = 0; i < numberGridCells; i++)
  {
    if (!neighborhoodList[i]) empty_cell++;
  }

  return ((double) (empty_cell) / (double) (numberGridCells));
}

Definition at line 584 of file PointGridSearch.cpp.

{
  // retrieve points over the current bounding box range

  facet_list.clean_out();
  DLIList<CubitPoint*> point_list;
  get_neighborhood_points( point_list );

  // retrieve all faces attached to the points in point_list

  for (int i = 0; i < point_list.size(); i++)
  {
    CubitPoint* point = point_list.get_and_step();

    DLIList<CubitFacet*> temp_facet_list;
    point->facets(temp_facet_list);

    for (int j = 0; j < temp_facet_list.size(); j++)
    {
      CubitFacet* facet = temp_facet_list.get_and_step();

      if (!facet->marked())
      {
        facet->marked(CUBIT_TRUE);
        facet_list.append(facet);
      }
    }
  }

  // unmark the found faces and return face_list

  for (int m = 0; m < facet_list.size(); m++)
  {
    facet_list.get_and_step()->marked(CUBIT_FALSE);
  }

}

Definition at line 501 of file PointGridSearch.cpp.

{
  // retrieve points over the current bounding box range

  point_list.clean_out();

  for (int k = boundingCellMinimumZ; k <= boundingCellMaximumZ; k++)
  {
    int kn = numberGridCellsY * k;
    for (int j = boundingCellMinimumY; j <= boundingCellMaximumY; j++)
    {
      int jn = numberGridCellsX * (kn + j);
      for (int i = boundingCellMinimumX; i <= boundingCellMaximumX; i++)
      {
        int in = jn + i;
        assert( in >= 0 && in < numberGridCells );
        if (neighborhoodList[in])
        {
           point_list += *(neighborhoodList[in]);
        }
      }
    }
  }
}
void PointGridSearch::get_neighborhood_points_sorted ( DLIList< CubitPoint * > &  point_list,
const CubitVector center,
double  cut_off 
)

Definition at line 535 of file PointGridSearch.cpp.

{
  point_list.clean_out();
  DLIList<CubitPoint*> temp_point_list;
  int i;

  for (int k = boundingCellMinimumZ; k <= boundingCellMaximumZ; k++)
  {
    int kn = numberGridCellsY * k;
    for (int j = boundingCellMinimumY; j <= boundingCellMaximumY; j++)
    {
      int jn = numberGridCellsX * (kn + j);
      for ( i = boundingCellMinimumX; i <= boundingCellMaximumX; i++)
      {
    int in = jn + i;
        if (neighborhoodList[in])
        {
          temp_point_list += *(neighborhoodList[in]);
        }
      }
    }
  }

  // evaluate point distance to center ... remove those larger than cut_off
  std::vector< std::pair<int, double> > sorted_index_list;
  CubitVector vec;
  cut_off *= cut_off;
  temp_point_list.reset();
  for ( i = 0; i < temp_point_list.size(); i++)
  {
    vec = center - temp_point_list.get_and_step()->coordinates();
    double distance = vec.length_squared();
    if (distance < cut_off)
    {
      sorted_index_list.push_back( std::make_pair( i, distance ) );
    }
  }

  std::sort(sorted_index_list.begin(), sorted_index_list.end(), sort_by_second<std::pair<int, double> >());
  temp_point_list.reset();
  for ( size_t j = 0; j < sorted_index_list.size(); j++ )
  {
    int index = sorted_index_list[j].first;
    point_list.append( temp_point_list.next( index ) );
  }
}
int PointGridSearch::grid_cell_index ( CubitPoint point) [inline, protected]

Definition at line 212 of file PointGridSearch.hpp.

{
  CubitVector range_vec = point->coordinates();
  range_vec -= gridRangeMinimum;

  int i = (int) (range_vec.x() * gridCellWidthInverse.x());
  int j = (int) (range_vec.y() * gridCellWidthInverse.y());
  int k = (int) (range_vec.z() * gridCellWidthInverse.z());

  return numberGridCellsX * (numberGridCellsY * k + j) + i;
}

Definition at line 189 of file PointGridSearch.hpp.

{return gridCellWidth;}

Definition at line 188 of file PointGridSearch.hpp.

{return gridRangeMaximum;}

Definition at line 187 of file PointGridSearch.hpp.

{return gridRangeMinimum;}
int PointGridSearch::in_grid ( const CubitVector position)

Definition at line 270 of file PointGridSearch.cpp.

{
  CubitVector range_vec;
  int i, j, k, index;

  range_vec = position - gridRangeMinimum;
  i = (int) (range_vec.x() * gridCellWidthInverse.x());
  j = (int) (range_vec.y() * gridCellWidthInverse.y());
  k = (int) (range_vec.z() * gridCellWidthInverse.z());
  index = numberGridCellsX * (numberGridCellsY * k + j) + i;
  if(index >= 0 && index < numberGridCells)
  {
    return CUBIT_TRUE;
  }
  else
  {
    return CUBIT_FALSE;
  }
}

Definition at line 193 of file PointGridSearch.hpp.

{return numberGridCells;}

Definition at line 190 of file PointGridSearch.hpp.

{return numberGridCellsX;}

Definition at line 191 of file PointGridSearch.hpp.

{return numberGridCellsY;}

Definition at line 192 of file PointGridSearch.hpp.

{return numberGridCellsZ;}

Definition at line 326 of file PointGridSearch.cpp.

{
    // remove point from the appropriate list
    int index = -1;
    ToolData* td_temp = point->get_TD(&TDCellIndex::is_cell_index);
    if (!td_temp)
      return;
    TDCellIndex* td_index = CAST_TO(td_temp, TDCellIndex);
    if( td_index ){
      index = td_index->cell_index();
      point->delete_TD( &TDCellIndex::is_cell_index );
    }
    remove_point_from_cell(point, index);
}
void PointGridSearch::remove_point_from_cell ( CubitPoint point,
int  i 
) [private]

Definition at line 679 of file PointGridSearch.cpp.

{
  if( i < 0 ) return;
  if (neighborhoodList[i])
  { 
    if (neighborhoodList[i]->move_to(point))
    {
      neighborhoodList[i]->extract();
    }
    if (!neighborhoodList[i]->size())
    {
      delete neighborhoodList[i];
      neighborhoodList[i] = NULL;
    }
  }
}

Definition at line 341 of file PointGridSearch.cpp.

{
    // set the range and cell min and max values for the vector

    boundingRangeMinimum = vec;
    boundingRangeMaximum = vec;

    cell_from_range();
}
void PointGridSearch::set_neighborhood_bounds ( const CubitVector center,
double  size 
)

Definition at line 351 of file PointGridSearch.cpp.

{
    // set the range and cell min and max values as a cube of breadth =
    // 2 * size centered at center

    CubitVector search_extent_box(size, size, size);
    boundingRangeMinimum = center - search_extent_box;
    boundingRangeMaximum = center + search_extent_box;

    cell_from_range();
}
void PointGridSearch::set_neighborhood_bounds ( CubitVector center,
double  x,
double  y,
double  z 
)

Definition at line 378 of file PointGridSearch.cpp.

{
    // set the range and cell min and max values as a box of size 2x, 2y, and
    // 2z, centered at center

    CubitVector search_extent_box(x, y, z);
    boundingRangeMinimum = center - search_extent_box;
    boundingRangeMaximum = center + search_extent_box;

    cell_from_range();
}

Definition at line 391 of file PointGridSearch.cpp.

{
    // initialize min and max range values to the first vectors coordinates
   
    boundingRangeMinimum = vec_1;
    boundingRangeMaximum = vec_1;

    // find the range and cell min and max values for the two vectors

    bounding_range(vec_2);
    cell_from_range();
}

Definition at line 405 of file PointGridSearch.cpp.

{
    // initialize min and max range values to the point coordinates

    boundingRangeMinimum = point->coordinates();
    boundingRangeMaximum = point->coordinates();

    // find the range and cell min and max values for the point and vector

    bounding_range(vec);
    cell_from_range();
}

Definition at line 418 of file PointGridSearch.cpp.

{
    // initialize min and max range values to the point coordinates

    boundingRangeMinimum = point->coordinates();
    boundingRangeMaximum = point->coordinates();

    // find the cell min and max values

    cell_from_range();
}
void PointGridSearch::set_neighborhood_bounds ( CubitPoint point_1,
CubitPoint point_2 
)

Definition at line 430 of file PointGridSearch.cpp.

{
    // initialize min and max range values to the first point coordinates

    boundingRangeMinimum = point_1->coordinates();
    boundingRangeMaximum = point_1->coordinates();

    // find the range and cell min and max values for the two points

    bounding_range(point_2);
    cell_from_range();
}

Definition at line 443 of file PointGridSearch.cpp.

{
    // find the range and cell min and max values for the point list

    bounding_range(point_list);
    cell_from_range();
}

Definition at line 451 of file PointGridSearch.cpp.

{
    // find the points belonging to the face

    DLIList<CubitPoint*> point_list;
    facet->points(point_list);

    // find the bounding range and cell min and max values for the point list

    bounding_range(point_list);
    cell_from_range();
}

Definition at line 464 of file PointGridSearch.cpp.

{
    DLIList<CubitPoint*> point_list;

    // find all points belonging to the faces in face_list
    int i;
    for ( i = 0; i < facet_list.size(); i++)
    {
      CubitFacet* facet = facet_list.get_and_step();

      DLIList<CubitPoint*> temp_point_list;
      facet->points(temp_point_list);
      for (int j = 0; j < temp_point_list.size(); j++)
      {
    CubitPoint* point = temp_point_list.get_and_step();

    if (!point->marked())
    {
      point->marked(CUBIT_TRUE);
      point_list.append(point);
    }
      }
    }

    // unmark the found points

    for (i = 0; i < point_list.size(); i++)
    {
      point_list.get_and_step()->marked(CUBIT_FALSE);
    }

    // find the bounding range and cell min and max values for the point list

    bounding_range(point_list);
    cell_from_range();
}
void PointGridSearch::set_neighborhood_bounds_max_edge ( const CubitVector center,
double  factor = 1.0 
)

Definition at line 363 of file PointGridSearch.cpp.

{
    // set the range and cell min and max values as a cube of breadth =
    // 2 * size centered at center

    CubitVector search_extent_box(maxEdgeLength*factor, 
                                  maxEdgeLength*factor, 
                                  maxEdgeLength*factor);
    boundingRangeMinimum = center - search_extent_box;
    boundingRangeMaximum = center + search_extent_box;

    cell_from_range();
}

Definition at line 651 of file PointGridSearch.cpp.

{
  int total_points = 0;

  for (int i = 0; i < numberGridCells; i++)
  {
    if (neighborhoodList[i])
    {
      total_points += neighborhoodList[i]->size();
    }
  }

  return total_points;
}

Member Data Documentation

Definition at line 78 of file PointGridSearch.hpp.

Definition at line 79 of file PointGridSearch.hpp.

Definition at line 80 of file PointGridSearch.hpp.

Definition at line 75 of file PointGridSearch.hpp.

Definition at line 76 of file PointGridSearch.hpp.

Definition at line 77 of file PointGridSearch.hpp.

Definition at line 64 of file PointGridSearch.hpp.

Definition at line 63 of file PointGridSearch.hpp.

Definition at line 62 of file PointGridSearch.hpp.

Definition at line 52 of file PointGridSearch.hpp.

Definition at line 49 of file PointGridSearch.hpp.

Definition at line 69 of file PointGridSearch.hpp.

Definition at line 66 of file PointGridSearch.hpp.

Definition at line 67 of file PointGridSearch.hpp.

Definition at line 68 of file PointGridSearch.hpp.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines