OGS
GeoLib::Grid< POINT > Class Template Referencefinal

Detailed Description

template<typename POINT>
class GeoLib::Grid< POINT >

Definition at line 21 of file Grid.h.

#include <Grid.h>

Inheritance diagram for GeoLib::Grid< POINT >:
[legend]
Collaboration diagram for GeoLib::Grid< POINT >:
[legend]

Public Member Functions

template<typename InputIterator>
 Grid (InputIterator first, InputIterator last, std::size_t max_num_per_grid_cell=1)
 The constructor of the grid object takes a vector of points or nodes. Furthermore the user can specify the average maximum number of points per grid cell.
 Grid (Grid const &)=delete
Gridoperator= (Grid const &)=delete
 ~Grid ()
template<typename P>
POINTgetNearestPoint (P const &pnt) const
template<typename P>
std::vector< std::size_t > getPointsInEpsilonEnvironment (P const &pnt, double eps) const
template<typename P>
std::vector< std::vector< POINT * > const * > getPntVecsOfGridCellsIntersectingCube (P const &center, double half_len) const
std::vector< std::vector< POINT * > const * > getPntVecsOfGridCellsIntersectingCuboid (Eigen::Vector3d const &min_pnt, Eigen::Vector3d const &max_pnt) const
void createGridGeometry (GeoLib::GEOObjects *geo_obj) const
Public Member Functions inherited from GeoLib::AABB
template<typename PNT_TYPE>
 AABB (std::vector< PNT_TYPE * > const &pnts, std::vector< std::size_t > const &ids)
template<typename InputIterator>
 AABB (InputIterator first, InputIterator last)
template<typename PNT_TYPE>
bool update (PNT_TYPE const &p)
template<typename T>
bool containsPoint (T const &pnt, double eps) const
template<typename T>
bool containsPointXY (T const &pnt) const
MinMaxPoints getMinMaxPoints () const
Eigen::Vector3d const & getMinPoint () const
Eigen::Vector3d const & getMaxPoint () const
bool containsAABB (AABB const &other_aabb) const

Private Member Functions

void initNumberOfSteps (std::size_t n_per_cell, std::size_t n_pnts, std::array< double, 3 > const &extensions)
template<typename T>
std::array< std::size_t, 3 > getGridCoords (T const &pnt) const
template<typename P>
std::array< double, 6 > getPointCellBorderDistances (P const &pnt, std::array< std::size_t, 3 > const &coords) const
template<typename P>
bool calcNearestPointInGridCell (P const &pnt, std::array< std::size_t, 3 > const &coords, double &sqr_min_dist, POINT *&nearest_pnt) const

Static Private Member Functions

static POINTcopyOrAddress (POINT &p)
static POINT const * copyOrAddress (POINT const &p)
static POINTcopyOrAddress (POINT *p)

Private Attributes

std::array< std::size_t, 3 > _n_steps = {{1, 1, 1}}
std::array< double, 3 > _step_sizes = {{0.0, 0.0, 0.0}}
std::vector< POINT * > * _grid_cell_nodes_map = nullptr

Constructor & Destructor Documentation

◆ Grid() [1/2]

template<typename POINT>
template<typename InputIterator>
GeoLib::Grid< POINT >::Grid ( InputIterator first,
InputIterator last,
std::size_t max_num_per_grid_cell = 1 )

The constructor of the grid object takes a vector of points or nodes. Furthermore the user can specify the average maximum number of points per grid cell.

The number of grid cells are computed with the following formula \(\frac{n_{points}}{n_{cells}} \le n_{max\_per\_cell}\)

In order to limit the memory wasting the maximum number of points per grid cell (in the average) should be a power of two (since std::vector objects resize itself with this step size).

Parameters
first,lastthe range of elements to examine
max_num_per_grid_cell(input) max number per grid cell in the average

Definition at line 184 of file Grid.h.

187{
188 auto const n_pnts(std::distance(first, last));
189
190 auto const& max{getMaxPoint()};
191 auto const& min{getMinPoint()};
193 {max[0] - min[0], max[1] - min[1], max[2] - min[2]}};
194
195 // enlarge delta
196 constexpr double direction = std::numeric_limits<double>::max();
198 [](double const d) { return std::nextafter(d, direction); });
199
200 assert(n_pnts > 0);
202 delta);
203
204 const std::size_t n_plane(_n_steps[0] * _n_steps[1]);
206
207 // some frequently used expressions to fill the grid vectors
208 for (std::size_t k(0); k < 3; k++)
209 {
211 {
213 }
214 _step_sizes[k] = delta[k] / _n_steps[k];
215 }
216
217 // fill the grid vectors
219 while (it != last)
220 {
222 if (coords < _n_steps)
223 {
224 std::size_t const pos(coords[0] + coords[1] * _n_steps[0] +
225 coords[2] * n_plane);
226 _grid_cell_nodes_map[pos].push_back(
227 const_cast<POINT*>(copyOrAddress(*it)));
228 }
229 else
230 {
231 ERR("Grid constructor: error computing indices [{:d}, {:d}, {:d}], "
232 "max indices [{:d}, {:d}, {:d}].",
233 coords[0], coords[1], coords[2], _n_steps[0], _n_steps[1],
234 _n_steps[2]);
235 }
236 it++;
237 }
238}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
Eigen::Vector3d const & getMaxPoint() const
Definition AABB.h:176
Eigen::Vector3d const & getMinPoint() const
Definition AABB.h:169
std::array< std::size_t, 3 > getGridCoords(T const &pnt) const
Definition Grid.h:408
static POINT * copyOrAddress(POINT &p)
Definition Grid.h:170
std::array< std::size_t, 3 > _n_steps
Definition Grid.h:174
std::array< double, 3 > _step_sizes
Definition Grid.h:175
std::vector< POINT * > * _grid_cell_nodes_map
Definition Grid.h:179
void initNumberOfSteps(std::size_t n_per_cell, std::size_t n_pnts, std::array< double, 3 > const &extensions)
Definition Grid.h:555

References GeoLib::AABB::AABB(), _grid_cell_nodes_map, _n_steps, _step_sizes, copyOrAddress(), ERR(), getGridCoords(), GeoLib::AABB::getMaxPoint(), GeoLib::AABB::getMinPoint(), initNumberOfSteps(), and GeoLib::POINT.

Referenced by Grid(), and operator=().

◆ Grid() [2/2]

template<typename POINT>
GeoLib::Grid< POINT >::Grid ( Grid< POINT > const & )
delete

References Grid().

◆ ~Grid()

template<typename POINT>
GeoLib::Grid< POINT >::~Grid ( )
inline

This is the destructor of the class. It deletes the internal data structures not including the pointers to the points.

Definition at line 51 of file Grid.h.

51{ delete[] _grid_cell_nodes_map; }

References _grid_cell_nodes_map.

Member Function Documentation

◆ calcNearestPointInGridCell()

template<typename POINT>
template<typename P>
bool GeoLib::Grid< POINT >::calcNearestPointInGridCell ( P const & pnt,
std::array< std::size_t, 3 > const & coords,
double & sqr_min_dist,
POINT *& nearest_pnt ) const
private

Definition at line 632 of file Grid.h.

637{
638 const std::size_t grid_idx(coords[0] + coords[1] * _n_steps[0] +
639 coords[2] * _n_steps[0] * _n_steps[1]);
643 if (pnts.empty())
644 return false;
645
646 const std::size_t n_pnts(pnts.size());
648 (pnts[0]->asEigenVector3d() - pnt.asEigenVector3d()).squaredNorm();
649 nearest_pnt = pnts[0];
650 for (std::size_t i(1); i < n_pnts; i++)
651 {
652 const double sqr_dist(
653 (pnts[i]->asEigenVector3d() - pnt.asEigenVector3d()).squaredNorm());
655 {
657 nearest_pnt = pnts[i];
658 }
659 }
660 return true;
661}

References _grid_cell_nodes_map, _n_steps, and GeoLib::POINT.

Referenced by getNearestPoint().

◆ copyOrAddress() [1/3]

template<typename POINT>
POINT * GeoLib::Grid< POINT >::copyOrAddress ( POINT & p)
inlinestaticprivate

Definition at line 170 of file Grid.h.

170{ return &p; }

References GeoLib::POINT.

Referenced by Grid().

◆ copyOrAddress() [2/3]

template<typename POINT>
POINT * GeoLib::Grid< POINT >::copyOrAddress ( POINT * p)
inlinestaticprivate

Definition at line 172 of file Grid.h.

172{ return p; }

References GeoLib::POINT.

◆ copyOrAddress() [3/3]

template<typename POINT>
POINT const * GeoLib::Grid< POINT >::copyOrAddress ( POINT const & p)
inlinestaticprivate

Definition at line 171 of file Grid.h.

171{ return &p; }

References GeoLib::POINT.

◆ createGridGeometry()

template<typename POINT>
void GeoLib::Grid< POINT >::createGridGeometry ( GeoLib::GEOObjects * geo_obj) const

Method creates a geometry for every mesh grid box. Additionally it creates one geometry containing all the box geometries.

Parameters
geo_obj

Definition at line 306 of file Grid.h.

307{
309
310 GeoLib::Point const& llf(getMinPoint());
311 GeoLib::Point const& urb(getMaxPoint());
312
313 const double dx((urb[0] - llf[0]) / _n_steps[0]);
314 const double dy((urb[1] - llf[1]) / _n_steps[1]);
315 const double dz((urb[2] - llf[2]) / _n_steps[2]);
316
317 // create grid names and grid boxes as geometry
318 for (std::size_t i(0); i < _n_steps[0]; i++)
319 {
320 for (std::size_t j(0); j < _n_steps[1]; j++)
321 {
322 for (std::size_t k(0); k < _n_steps[2]; k++)
323 {
324 std::string name("Grid-");
326 name += "-";
328 name += "-";
330 grid_names.push_back(name);
331
332 {
334 points.push_back(new GeoLib::Point(
335 llf[0] + i * dx, llf[1] + j * dy, llf[2] + k * dz));
336 points.push_back(new GeoLib::Point(llf[0] + i * dx,
337 llf[1] + (j + 1) * dy,
338 llf[2] + k * dz));
339 points.push_back(new GeoLib::Point(llf[0] + (i + 1) * dx,
340 llf[1] + (j + 1) * dy,
341 llf[2] + k * dz));
342 points.push_back(new GeoLib::Point(llf[0] + (i + 1) * dx,
343 llf[1] + j * dy,
344 llf[2] + k * dz));
345 points.push_back(new GeoLib::Point(llf[0] + i * dx,
346 llf[1] + j * dy,
347 llf[2] + (k + 1) * dz));
348 points.push_back(new GeoLib::Point(llf[0] + i * dx,
349 llf[1] + (j + 1) * dy,
350 llf[2] + (k + 1) * dz));
351 points.push_back(new GeoLib::Point(llf[0] + (i + 1) * dx,
352 llf[1] + (j + 1) * dy,
353 llf[2] + (k + 1) * dz));
354 points.push_back(new GeoLib::Point(llf[0] + (i + 1) * dx,
355 llf[1] + j * dy,
356 llf[2] + (k + 1) * dz));
357 geo_obj->addPointVec(std::move(points), grid_names.back());
358 }
359
361 auto const& points = *geo_obj->getPointVec(grid_names.back());
362
363 auto* ply0(new GeoLib::Polyline(points));
364 for (std::size_t l(0); l < 4; l++)
365 ply0->addPoint(l);
366 ply0->addPoint(0);
367 plys.push_back(ply0);
368
369 auto* ply1(new GeoLib::Polyline(points));
370 for (std::size_t l(4); l < 8; l++)
371 ply1->addPoint(l);
372 ply1->addPoint(4);
373 plys.push_back(ply1);
374
375 auto* ply2(new GeoLib::Polyline(points));
376 ply2->addPoint(0);
377 ply2->addPoint(4);
378 plys.push_back(ply2);
379
380 auto* ply3(new GeoLib::Polyline(points));
381 ply3->addPoint(1);
382 ply3->addPoint(5);
383 plys.push_back(ply3);
384
385 auto* ply4(new GeoLib::Polyline(points));
386 ply4->addPoint(2);
387 ply4->addPoint(6);
388 plys.push_back(ply4);
389
390 auto* ply5(new GeoLib::Polyline(points));
391 ply5->addPoint(3);
392 ply5->addPoint(7);
393 plys.push_back(ply5);
394
395 geo_obj->addPolylineVec(std::move(plys), grid_names.back(),
397 }
398 }
399 }
401
402 geo_obj->mergeGeometries(grid_names, merged_geo_name);
403}
std::map< std::string, std::size_t > NameIdMap
Definition TemplateVec.h:30

References _n_steps, GeoLib::GEOObjects::addPointVec(), GeoLib::GEOObjects::addPolylineVec(), GeoLib::AABB::getMaxPoint(), GeoLib::AABB::getMinPoint(), GeoLib::GEOObjects::getPointVec(), and GeoLib::GEOObjects::mergeGeometries().

◆ getGridCoords()

template<typename POINT>
template<typename T>
std::array< std::size_t, 3 > GeoLib::Grid< POINT >::getGridCoords ( T const & pnt) const
private

Method calculates the grid cell coordinates for the given point pnt. If the point is located outside of the bounding box the coordinates of the grid cell on the border that is nearest to given point will be returned.

Parameters
pnt(input) the coordinates of the point
Returns
the coordinates of the grid cell

Definition at line 408 of file Grid.h.

409{
410 auto const [min_point, max_point] = getMinMaxPoints();
411
413 for (std::size_t k(0); k < 3; k++)
414 {
415 if (pnt[k] < min_point[k])
416 {
417 continue;
418 }
419 if (pnt[k] >= max_point[k])
420 {
421 coords[k] = _n_steps[k] - 1;
422 continue;
423 }
424 coords[k] = static_cast<std::size_t>(
425 std::floor((pnt[k] - min_point[k]) /
428 }
429 return coords;
430}
MinMaxPoints getMinMaxPoints() const
Definition AABB.h:163

References _n_steps, _step_sizes, and GeoLib::AABB::getMinMaxPoints().

Referenced by Grid(), getNearestPoint(), getPntVecsOfGridCellsIntersectingCube(), and getPntVecsOfGridCellsIntersectingCuboid().

◆ getNearestPoint()

template<typename POINT>
template<typename P>
POINT * GeoLib::Grid< POINT >::getNearestPoint ( P const & pnt) const

The method calculates the grid cell the given point is belonging to, i.e., the (internal) coordinates of the grid cell are computed. The method searches the actual grid cell and all its neighbors for the POINT object which has the smallest distance. A pointer to this object is returned.

If there is not such a point, i.e., all the searched grid cells do not contain any POINT object a nullptr is returned.

Parameters
pntsearch point
Returns
a pointer to the point with the smallest distance within the grid cells that are outlined above or nullptr

Definition at line 458 of file Grid.h.

459{
461 auto const& min_point{getMinPoint()};
462 auto const& max_point{getMaxPoint()};
463
465 POINT* nearest_pnt(nullptr);
466
468
470 {
472 if (dists[0] >= min_dist && dists[1] >= min_dist &&
473 dists[2] >= min_dist && dists[3] >= min_dist &&
474 dists[4] >= min_dist && dists[5] >= min_dist)
475 {
476 return nearest_pnt;
477 }
478 }
479 else
480 {
481 // search in all border cells for at least one neighbor
482 double sqr_min_dist_tmp;
483 POINT* nearest_pnt_tmp(nullptr);
485
486 while (nearest_pnt == nullptr)
487 {
489 {coords[0] < offset ? 0 : coords[0] - offset,
490 coords[1] < offset ? 0 : coords[1] - offset,
491 coords[2] < offset ? 0 : coords[2] - offset}};
492 for (; ijk[0] < coords[0] + offset; ijk[0]++)
493 {
494 for (; ijk[1] < coords[1] + offset; ijk[1]++)
495 {
496 for (; ijk[2] < coords[2] + offset; ijk[2]++)
497 {
498 // do not check the origin grid cell twice
499 if (ijk[0] == coords[0] && ijk[1] == coords[1] &&
500 ijk[2] == coords[2])
501 {
502 continue;
503 }
504 // check if temporary grid cell coordinates are valid
505 if (ijk[0] >= _n_steps[0] || ijk[1] >= _n_steps[1] ||
506 ijk[2] >= _n_steps[2])
507 {
508 continue;
509 }
510
513 {
515 {
518 }
519 }
520 } // end k
521 } // end j
522 } // end i
523 offset++;
524 } // end while
525 } // end else
526
527 double const len(
528 (pnt.asEigenVector3d() - nearest_pnt->asEigenVector3d()).norm());
529 // search all other grid cells within the cube with the edge nodes
532
533 const std::size_t n_vecs(vecs_of_pnts.size());
534 for (std::size_t j(0); j < n_vecs; j++)
535 {
537 const std::size_t n_pnts(pnts.size());
538 for (std::size_t k(0); k < n_pnts; k++)
539 {
540 const double sqr_dist(
541 (pnt.asEigenVector3d() - pnts[k]->asEigenVector3d())
542 .squaredNorm());
544 {
546 nearest_pnt = pnts[k];
547 }
548 }
549 }
550
551 return nearest_pnt;
552}
std::array< double, 6 > getPointCellBorderDistances(P const &pnt, std::array< std::size_t, 3 > const &coords) const
Definition Grid.h:434
std::vector< std::vector< POINT * > const * > getPntVecsOfGridCellsIntersectingCube(P const &center, double half_len) const
Definition Grid.h:243
bool calcNearestPointInGridCell(P const &pnt, std::array< std::size_t, 3 > const &coords, double &sqr_min_dist, POINT *&nearest_pnt) const
Definition Grid.h:632

References _n_steps, calcNearestPointInGridCell(), getGridCoords(), GeoLib::AABB::getMaxPoint(), GeoLib::AABB::getMinPoint(), getPntVecsOfGridCellsIntersectingCube(), getPointCellBorderDistances(), and GeoLib::POINT.

◆ getPntVecsOfGridCellsIntersectingCube()

template<typename POINT>
template<typename P>
std::vector< std::vector< POINT * > const * > GeoLib::Grid< POINT >::getPntVecsOfGridCellsIntersectingCube ( P const & center,
double half_len ) const

Method fetches the vectors of all grid cells intersecting the axis aligned cuboid defined by two points. The first point with minimal coordinates in all directions. The second point with maximal coordinates in all directions.

Parameters
center(input) the center point of the axis aligned cube
half_len(input) half of the edge length of the axis aligned cube
Returns
vector of vectors of points within grid cells that intersects the axis aligned cube

Definition at line 243 of file Grid.h.

245{
246 Eigen::Vector3d const c{center[0], center[1], center[2]};
249
251 pnts.reserve(
254 .prod());
255
257 for (std::size_t c0 = min_coords[0]; c0 < max_coords[0] + 1; c0++)
258 {
259 for (std::size_t c1 = min_coords[1]; c1 < max_coords[1] + 1; c1++)
260 {
262 for (std::size_t c2 = min_coords[2]; c2 < max_coords[2] + 1; c2++)
263 {
264 pnts.push_back(
266 c2 * steps0_x_steps1]));
267 }
268 }
269 }
270 return pnts;
271}

References _grid_cell_nodes_map, _n_steps, and getGridCoords().

Referenced by MeshToolsLib::MeshRevision::collapseNodeIndices(), getNearestPoint(), and getPointsInEpsilonEnvironment().

◆ getPntVecsOfGridCellsIntersectingCuboid()

template<typename POINT>
std::vector< std::vector< POINT * > const * > GeoLib::Grid< POINT >::getPntVecsOfGridCellsIntersectingCuboid ( Eigen::Vector3d const & min_pnt,
Eigen::Vector3d const & max_pnt ) const

◆ getPointCellBorderDistances()

template<typename POINT>
template<typename P>
std::array< double, 6 > GeoLib::Grid< POINT >::getPointCellBorderDistances ( P const & pnt,
std::array< std::size_t, 3 > const & coords ) const
private

point numbering of the grid cell is as follow

7 -------- 6
/: /|
/ : / |
/ : / |
/ : / |
4 -------- 5 |
| 3 ....|... 2
| . | /
| . | /
| . | /
|. |/
0 -------- 1

the face numbering is as follow: face nodes 0 0,3,2,1 bottom 1 0,1,5,4 front 2 1,2,6,5 right 3 2,3,7,6 back 4 3,0,4,7 left 5 4,5,6,7 top

Parameters
pnt(input) coordinates of the point
coordsof the grid cell
Returns
squared distances of the point to the faces of the grid cell ordered in the same sequence as above described

Definition at line 434 of file Grid.h.

436{
437 auto const& min_point{getMinPoint()};
439 dists[0] =
440 std::abs(pnt[2] - min_point[2] + coords[2] * _step_sizes[2]); // bottom
441 dists[5] = std::abs(pnt[2] - min_point[2] +
442 (coords[2] + 1) * _step_sizes[2]); // top
443
444 dists[1] =
445 std::abs(pnt[1] - min_point[1] + coords[1] * _step_sizes[1]); // front
446 dists[3] = std::abs(pnt[1] - min_point[1] +
447 (coords[1] + 1) * _step_sizes[1]); // back
448
449 dists[4] =
450 std::abs(pnt[0] - min_point[0] + coords[0] * _step_sizes[0]); // left
451 dists[2] = std::abs(pnt[0] - min_point[0] +
452 (coords[0] + 1) * _step_sizes[0]); // right
453 return dists;
454}

References _step_sizes, and GeoLib::AABB::getMinPoint().

Referenced by getNearestPoint().

◆ getPointsInEpsilonEnvironment()

template<typename POINT>
template<typename P>
std::vector< std::size_t > GeoLib::Grid< POINT >::getPointsInEpsilonEnvironment ( P const & pnt,
double eps ) const

Definition at line 665 of file Grid.h.

667{
670
671 double const sqr_eps(eps * eps);
672
674 for (auto vec : vec_pnts)
675 {
676 for (auto const p : *vec)
677 {
678 if ((p->asEigenVector3d() - pnt.asEigenVector3d()).squaredNorm() <=
679 sqr_eps)
680 {
681 pnts.push_back(p->getID());
682 }
683 }
684 }
685
686 return pnts;
687}

References getPntVecsOfGridCellsIntersectingCube().

Referenced by MeshGeoToolsLib::MeshNodesOnPoint::MeshNodesOnPoint().

◆ initNumberOfSteps()

template<typename POINT>
void GeoLib::Grid< POINT >::initNumberOfSteps ( std::size_t n_per_cell,
std::size_t n_pnts,
std::array< double, 3 > const & extensions )
private

Computes the number of grid cells per spatial dimension the objects (points or mesh nodes) will be sorted in. On the one hand the number of grid cells should be small to reduce the management overhead. On the other hand the number should be large such that each grid cell contains only a small number of objects. Under the assumption that the points are distributed equidistant in space the grid cells should be as cubical as possible. At first it is necessary to determine the spatial dimension the grid should have. The dimensions are computed from the spatial extensions Let \(\max\) be the largest spatial extension. The grid will have a spatial dimension if the ratio of the corresponding spatial extension and the maximal extension is \(\ge 10^{-4}\). The second step consists of computing the number of cells per dimension.

Definition at line 555 of file Grid.h.

557{
558 double const max_extension(
559 *std::max_element(extensions.cbegin(), extensions.cend()));
560
561 std::bitset<3> dim; // all bits set to zero
562 for (std::size_t k(0); k < 3; ++k)
563 {
564 // set dimension if the ratio kth-extension/max_extension >= 1e-4
565 if (extensions[k] >= 1e-4 * max_extension)
566 {
567 dim[k] = true;
568 }
569 }
570
571 // structured grid: n_cells = _n_steps[0] * _n_steps[1] * _n_steps[2]
572 // *** condition: n_pnts / n_cells < n_per_cell
573 // => n_pnts / n_per_cell < n_cells
574 // _n_steps[1] = _n_steps[0] * extensions[1]/extensions[0],
575 // _n_steps[2] = _n_steps[0] * extensions[2]/extensions[0],
576 // => n_cells = _n_steps[0]^3 * extensions[1]/extensions[0] *
577 // extensions[2]/extensions[0],
578 // => _n_steps[0] = cbrt(n_cells * extensions[0]^2 /
579 // (extensions[1]*extensions[2]))
580 auto sc_ceil = [](double v) { return static_cast<std::size_t>(ceil(v)); };
581
582 switch (dim.count())
583 {
584 case 3: // 3d case
585 _n_steps[0] = sc_ceil(
587 (extensions[0] / extensions[2]) / n_per_cell));
588 _n_steps[1] = sc_ceil(
589 _n_steps[0] * std::min(extensions[1] / extensions[0], 100.0));
590 _n_steps[2] = sc_ceil(
591 _n_steps[0] * std::min(extensions[2] / extensions[0], 100.0));
592 break;
593 case 2: // 2d cases
594 if (dim[0] && dim[1])
595 { // xy
597 (n_per_cell * extensions[1])));
598 _n_steps[1] =
599 sc_ceil(_n_steps[0] *
600 std::min(extensions[1] / extensions[0], 100.0));
601 }
602 else if (dim[0] && dim[2])
603 { // xz
605 (n_per_cell * extensions[2])));
606 _n_steps[2] =
607 sc_ceil(_n_steps[0] *
608 std::min(extensions[2] / extensions[0], 100.0));
609 }
610 else if (dim[1] && dim[2])
611 { // yz
613 (n_per_cell * extensions[2])));
614 _n_steps[2] =
615 sc_ceil(std::min(extensions[2] / extensions[1], 100.0));
616 }
617 break;
618 case 1: // 1d cases
619 for (std::size_t k(0); k < 3; ++k)
620 {
621 if (dim[k])
622 {
623 _n_steps[k] =
624 sc_ceil(static_cast<double>(n_pnts) / n_per_cell);
625 }
626 }
627 }
628}

References _n_steps.

Referenced by Grid().

◆ operator=()

template<typename POINT>
Grid & GeoLib::Grid< POINT >::operator= ( Grid< POINT > const & )
delete

References Grid().

Member Data Documentation

◆ _grid_cell_nodes_map

template<typename POINT>
std::vector<POINT*>* GeoLib::Grid< POINT >::_grid_cell_nodes_map = nullptr
private

This is an array that stores pointers to POINT objects.

Definition at line 179 of file Grid.h.

Referenced by Grid(), ~Grid(), calcNearestPointInGridCell(), getPntVecsOfGridCellsIntersectingCube(), and getPntVecsOfGridCellsIntersectingCuboid().

◆ _n_steps

template<typename POINT>
std::array<std::size_t, 3> GeoLib::Grid< POINT >::_n_steps = {{1, 1, 1}}
private

◆ _step_sizes

template<typename POINT>
std::array<double, 3> GeoLib::Grid< POINT >::_step_sizes = {{0.0, 0.0, 0.0}}
private

Definition at line 175 of file Grid.h.

175{{0.0, 0.0, 0.0}};

Referenced by Grid(), getGridCoords(), and getPointCellBorderDistances().


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