Loading...
Searching...
No Matches
surface_skd_tree.cpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25#include "bitpit_common.hpp"
26
27#include "surface_skd_tree.hpp"
28
29namespace bitpit {
30
44SurfaceSkdTree::SurfaceSkdTree(const SurfaceKernel *patch, bool interiorCellsOnly)
45 : PatchSkdTree(patch, interiorCellsOnly)
46{
47}
48
55void SurfaceSkdTree::clear(bool release)
56{
57 if (release) {
58 m_closestCellCandidates = ClosestCellCandidates();
59 }
60
61 PatchSkdTree::clear(release);
62}
63
73double SurfaceSkdTree::evalPointDistance(const std::array<double, 3> &point) const
74{
75 return evalPointDistance(point, std::numeric_limits<double>::max(), false);
76}
77
94double SurfaceSkdTree::evalPointDistance(const std::array<double, 3> &point, double maxDistance) const
95{
96 return evalPointDistance(point, maxDistance, false);
97}
98
118double SurfaceSkdTree::evalPointDistance(const std::array<double, 3> &point, double maxDistance, bool interiorCellsOnly) const
119{
120 long id;
121 double distance = maxDistance;
122
123 findPointClosestCell(point, interiorCellsOnly, &id, &distance);
124
125 return distance;
126}
127
138void SurfaceSkdTree::evalPointDistance(int nPoints, const std::array<double, 3> *points, double *distances) const
139{
140 std::vector<double> maxDistances(nPoints, std::numeric_limits<double>::max());
141
142 evalPointDistance(nPoints, points, maxDistances.data(), false, distances);
143}
144
162void SurfaceSkdTree::evalPointDistance(int nPoints, const std::array<double, 3> *points, double maxDistance, double *distances) const
163{
164 std::vector<double> maxDistances(nPoints, maxDistance);
165
166 evalPointDistance(nPoints, points, maxDistances.data(), false, distances);
167}
168
186void SurfaceSkdTree::evalPointDistance(int nPoints, const std::array<double, 3> *points, const double *maxDistances, double *distances) const
187{
188 evalPointDistance(nPoints, points, maxDistances, false, distances);
189}
190
191
211void SurfaceSkdTree::evalPointDistance(int nPoints, const std::array<double, 3> *points, const double *maxDistances, bool interiorCellsOnly, double *distances) const
212{
213 std::vector<long> ids(nPoints, Cell::NULL_ID);
214
215 findPointClosestCell(nPoints, points, maxDistances, interiorCellsOnly, ids.data(), distances);
216}
217
218#if BITPIT_ENABLE_MPI
230void SurfaceSkdTree::evalPointGlobalDistance(int nPoints, const std::array<double, 3> *points, double *distances) const
231{
232 std::vector<double> maxDistances(nPoints, std::numeric_limits<double>::max());
233
234 evalPointGlobalDistance(nPoints, points, maxDistances.data(), distances);
235}
236
255void SurfaceSkdTree::evalPointGlobalDistance(int nPoints, const std::array<double, 3> *points, double maxDistance, double *distances) const
256{
257 std::vector<double> maxDistances(nPoints, maxDistance);
258
259 evalPointGlobalDistance(nPoints, points, maxDistances.data(), distances);
260}
261
280void SurfaceSkdTree::evalPointGlobalDistance(int nPoints, const std::array<double, 3> *points, const double *maxDistances, double *distances) const
281{
282 std::vector<long> ids(nPoints, Cell::NULL_ID);
283 std::vector<int> ranks(nPoints, -1);
284
285 findPointClosestGlobalCell(nPoints, points, maxDistances, ids.data(), ranks.data(), distances);
286}
287#endif
288
304long SurfaceSkdTree::findPointClosestCell(const std::array<double, 3> &point, long *id, double *distance) const
305{
306 return findPointClosestCell(point, std::numeric_limits<double>::max(), false, id, distance);
307}
308
326long SurfaceSkdTree::findPointClosestCell(const std::array<double, 3> &point, double maxDistance,
327 long *id, double *distance) const
328{
329 return findPointClosestCell(point, maxDistance, false, id, distance);
330}
331
351long SurfaceSkdTree::findPointClosestCell(const std::array<double, 3> &point, double maxDistance,
352 bool interiorCellsOnly, long *id, double *distance) const
353{
354 // Tolerance for distance evaluations
355 const PatchKernel &patch = getPatch();
356 double tolerance = patch.getTol();
357
358 // Initialize the cell id
359 *id = Cell::NULL_ID;
360
361 // Get the root of the tree
362 std::size_t rootId = 0;
363 const SkdNode &root = m_nodes[rootId];
364 if (root.isEmpty()) {
365 *distance = std::numeric_limits<double>::max();
366
367 return 0;
368 }
369
370 ClosestCellCandidates *closestCellCandidates = nullptr;
371 std::unique_ptr<ClosestCellCandidates> privateStorage = nullptr;
372 if (areLookupsThreadSafe()) {
373 privateStorage = std::unique_ptr<ClosestCellCandidates>(new ClosestCellCandidates());
374 closestCellCandidates = privateStorage.get();
375 } else {
376 closestCellCandidates = &m_closestCellCandidates;
377 }
378 closestCellCandidates->maxDistance = *distance;
379 findPointClosestCandidates(point, maxDistance, closestCellCandidates);
380
381 // Process the candidates and find the closest cell
382 long nDistanceEvaluations = 0;
383 for (std::size_t k = 0; k < closestCellCandidates->ids.size(); ++k) {
384 // Do not consider nodes with a minimum distance greater than the
385 // distance estimate
386 if (utils::DoubleFloatingGreater()(closestCellCandidates->minDistances.at(k), closestCellCandidates->maxDistance, tolerance, tolerance)) {
387 continue;
388 }
389
390 // Evaluate the distance
391 std::size_t nodeId = closestCellCandidates->ids.at(k);
392 const SkdNode &node = m_nodes[nodeId];
393
394 node.updatePointClosestCell(point, interiorCellsOnly, id, &(closestCellCandidates->maxDistance));
395 ++nDistanceEvaluations;
396 }
397 *distance = closestCellCandidates->maxDistance;
398
399 return nDistanceEvaluations;
400}
401
413long SurfaceSkdTree::findPointClosestCell(int nPoints, const std::array<double, 3> *points, long *ids, double *distances) const
414{
415 long nDistanceEvaluations = 0;
416 for (int i = 0; i < nPoints; ++i) {
417 nDistanceEvaluations += findPointClosestCell(points[i], std::numeric_limits<double>::max(), false, ids + i, distances + i);
418 }
419
420 return nDistanceEvaluations;
421}
422
439long SurfaceSkdTree::findPointClosestCell(int nPoints, const std::array<double, 3> *points, double maxDistance, long *ids, double *distances) const
440{
441 long nDistanceEvaluations = 0;
442 for (int i = 0; i < nPoints; ++i) {
443 nDistanceEvaluations += findPointClosestCell(points[i], maxDistance, false, ids + i, distances + i);
444 }
445
446 return nDistanceEvaluations;
447}
448
466long SurfaceSkdTree::findPointClosestCell(int nPoints, const std::array<double, 3> *points, const double *maxDistances, long *ids, double *distances) const
467{
468 return findPointClosestCell(nPoints, points, maxDistances, false, ids, distances);
469}
470
490long SurfaceSkdTree::findPointClosestCell(int nPoints, const std::array<double, 3> *points, const double *maxDistances, bool interiorCellsOnly, long *ids, double *distances) const
491{
492 long nDistanceEvaluations = 0;
493 for (int i = 0; i < nPoints; ++i) {
494 nDistanceEvaluations += findPointClosestCell(points[i], maxDistances[i], interiorCellsOnly, ids + i, distances + i);
495 }
496
497 return nDistanceEvaluations;
498}
499
500
515long SurfaceSkdTree::findPointClosestCells(const std::array<double, 3> &point, std::vector<long> &ids, double *distance) const
516{
517 return findPointClosestCells(point, std::numeric_limits<double>::max(), false, ids, distance);
518}
519
537long SurfaceSkdTree::findPointClosestCells(const std::array<double, 3> &point, double maxDistance,
538 std::vector<long> &ids, double *distance) const
539{
540 return findPointClosestCells(point, maxDistance, false, ids, distance);
541}
542
563long SurfaceSkdTree::findPointClosestCells(const std::array<double, 3> &point, double maxDistance,
564 bool interiorCellsOnly, std::vector<long> &ids, double *distance) const
565{
566 // Tolerance for distance evaluations
567 const PatchKernel &patch = getPatch();
568 double tolerance = patch.getTol();
569
570 // Initialize the cell id
571 ids.resize(1, Cell::NULL_ID);
572
573 // Get the root of the tree
574 std::size_t rootId = 0;
575 const SkdNode &root = m_nodes[rootId];
576 if (root.isEmpty()) {
577 *distance = std::numeric_limits<double>::max();
578
579 return 0;
580 }
581
582 // If threads safe lookups are not needed, a temporary data structure
583 // are declared as member of the class to avoid its reallocation every
584 // time the function is called.
585 ClosestCellCandidates *closestCellCandidates = nullptr;
586 std::unique_ptr<ClosestCellCandidates> privateStorage = nullptr;
587 if (areLookupsThreadSafe()) {
588 privateStorage = std::unique_ptr<ClosestCellCandidates>(new ClosestCellCandidates());
589 closestCellCandidates = privateStorage.get();
590 } else {
591 closestCellCandidates = &m_closestCellCandidates;
592 }
593 findPointClosestCandidates(point, maxDistance, closestCellCandidates);
594
595 // Process the candidates and find the closest cells
596 *distance = closestCellCandidates->maxDistance;
597 long nDistanceEvaluations = 0;
598 for (std::size_t k = 0; k < closestCellCandidates->ids.size(); ++k) {
599 // Do not consider nodes with a minimum distance greater than the
600 // distance estimate
601 if (utils::DoubleFloatingGreater()(closestCellCandidates->minDistances.at(k), *distance, tolerance, tolerance)) {
602 continue;
603 }
604
605 // Evaluate the distance
606 std::size_t nodeId = closestCellCandidates->ids.at(k);
607 const SkdNode &node = m_nodes[nodeId];
608
609 node.updatePointClosestCells(point, interiorCellsOnly, &ids, distance);
610 ++nDistanceEvaluations;
611 }
612
613 return nDistanceEvaluations;
614}
615
631void SurfaceSkdTree::findPointClosestCandidates(const std::array<double, 3> &point, double maxDistance,
632 ClosestCellCandidates *candidates) const
633{
634 // Clear candidates
635 candidates->clear();
636
637 // Tolerance for distance evaluations
638 const PatchKernel &patch = getPatch();
639 double tolerance = patch.getTol();
640
641 // Get the root of the tree
642 std::size_t rootId = 0;
643 const SkdNode &root = m_nodes[rootId];
644
645 // Initialize a distance estimate
646 //
647 // The real distance will be lesser than or equal to the estimate.
648 //
649 // Care must be taken to avoid overflow when performing the multiplication.
650 double squaredMaxDistance;
651 if (maxDistance <= 1. || maxDistance < std::numeric_limits<double>::max() / maxDistance) {
652 squaredMaxDistance = maxDistance * maxDistance;
653 } else {
654 squaredMaxDistance = std::numeric_limits<double>::max();
655 }
656
657 double squareDistanceEstimate = std::min(root.evalPointMaxSquareDistance(point), squaredMaxDistance);
658
659 // Get a list of candidates nodes
660 //
661 // First, we gather all the candidates and then we evaluate the distance
662 // of each candidate. Since distance estimate is constantly updated when
663 // new nodes are processed, the final estimate may be smaller than the
664 // minimum distance of some candidates. Processing the candidates after
665 // scanning all the tree, allows to discard some of them without the need
666 // of evaluating the exact distance.
667
668 candidates->nodeStack.push_back(rootId);
669 while (!candidates->nodeStack.empty()) {
670 std::size_t nodeId = candidates->nodeStack.back();
671 const SkdNode &node = m_nodes[nodeId];
672 candidates->nodeStack.pop_back();
673
674 // Do not consider nodes with a minimum distance greater than
675 // the distance estimate
676 double nodeMinSquareDistance = node.evalPointMinSquareDistance(point);
677 if (utils::DoubleFloatingGreater()(nodeMinSquareDistance, squareDistanceEstimate, tolerance, tolerance)) {
678 continue;
679 }
680
681 // Update the distance estimate
682 //
683 // The real distance will be less than or equal to the estimate.
684 double nodeMaxSquareDistance = node.evalPointMaxSquareDistance(point);
685 squareDistanceEstimate = std::min(nodeMaxSquareDistance, squareDistanceEstimate);
686
687 // If the node is a leaf add it to the candidates, otherwise add its
688 // children to the stack.
689 bool isLeaf = true;
690 for (int i = SkdNode::CHILD_BEGIN; i != SkdNode::CHILD_END; ++i) {
691 SkdNode::ChildLocation childLocation = static_cast<SkdNode::ChildLocation>(i);
692 std::size_t childId = node.getChildId(childLocation);
693 if (childId != SkdNode::NULL_ID) {
694 isLeaf = false;
695 candidates->nodeStack.push_back(childId);
696 }
697 }
698
699 if (isLeaf) {
700 candidates->ids.push_back(nodeId);
701 candidates->minDistances.push_back(std::sqrt(nodeMinSquareDistance));
702 }
703 }
704
705 // Set the distance estimate
706 if (!candidates->ids.empty()) {
707 candidates->maxDistance = std::sqrt(squareDistanceEstimate);
708 } else {
709 candidates->maxDistance = std::numeric_limits<double>::max();
710 }
711}
712
713#if BITPIT_ENABLE_MPI
728long SurfaceSkdTree::findPointClosestGlobalCell(int nPoints, const std::array<double, 3> *points,
729 long *ids, int *ranks, double *distances) const
730{
731 return findPointClosestGlobalCell(nPoints, points, std::numeric_limits<double>::max(), ids, ranks, distances);
732}
733
753long SurfaceSkdTree::findPointClosestGlobalCell(int nPoints, const std::array<double, 3> *points, double maxDistance,
754 long *ids, int *ranks, double *distances) const
755{
756 std::vector<double> maxDistances(nPoints, maxDistance);
757
758 return findPointClosestGlobalCell(nPoints, points, maxDistances.data(), ids, ranks, distances);
759}
760
781long SurfaceSkdTree::findPointClosestGlobalCell(int nPoints, const std::array<double, 3> *points, const double *maxDistances,
782 long *ids, int *ranks, double *distances) const
783{
784 long nDistanceEvaluations = 0;
785
786 // Early return is the patch is not partitioned
787 const PatchKernel &patch = getPatch();
788 if (!patch.isPartitioned()) {
789 for (int i = 0; i < nPoints; ++i) {
790 // Evaluate distance
791 nDistanceEvaluations += findPointClosestCell(points[i], maxDistances[i], ids + i, distances + i);
792
793 // The patch is not partitioned, all cells are local
794 ranks[i] = patch.getRank();
795 }
796
797 return nDistanceEvaluations;
798 }
799
800 // Get MPI communicator
801 if (!isCommunicatorSet()) {
802 throw std::runtime_error("Skd-tree communicator has not been set.");
803 }
804
805 MPI_Comm communicator = getCommunicator();
806
807 // Gather the number of points associated to each process
808 std::vector<int> pointsCount(m_nProcessors);
809 MPI_Allgather(&nPoints, 1, MPI_INT, pointsCount.data(), 1, MPI_INT, communicator);
810
811 // Evaluate information for data communications
812 std::vector<int> globalPointsDispls(m_nProcessors, 0);
813 std::vector<int> globalPointsOffsets(m_nProcessors, 0);
814 std::vector<int> globalPointsDataCount(m_nProcessors, 0);
815
816 globalPointsDataCount[0] = 3 * pointsCount[0];
817 for (int i = 1; i < m_nProcessors; ++i) {
818 globalPointsDispls[i] = globalPointsDispls[i - 1] + 3 * pointsCount[i - 1];
819 globalPointsOffsets[i] = globalPointsOffsets[i - 1] + pointsCount[i - 1];
820 globalPointsDataCount[i] = 3 * pointsCount[i];
821 }
822
823 int nGlobalPoints = globalPointsDispls.back() + pointsCount.back();
824
825 // Gather point coordinates
826 std::vector<std::array<double,3>> globalPoints(nGlobalPoints);
827 int pointsDataCount = 3 * nPoints;
828 MPI_Allgatherv(points, pointsDataCount, MPI_DOUBLE, globalPoints.data(),
829 globalPointsDataCount.data(), globalPointsDispls.data(), MPI_DOUBLE, communicator);
830
831 // Gather vector with all maximum distances
832 std::vector<double> globalMaxDistances(nGlobalPoints);
833 MPI_Allgatherv(maxDistances, nPoints, MPI_DOUBLE, globalMaxDistances.data(),
834 pointsCount.data(), globalPointsOffsets.data(), MPI_DOUBLE, communicator);
835
836 // Initialize distance information
837 std::vector<SkdGlobalCellDistance> globalCellDistances(nGlobalPoints);
838
839 // Call local find point closest cell for each global point collected
840 for (int i = 0; i < nGlobalPoints; ++i) {
841 // Get point information
842 const std::array<double, 3> &point = globalPoints[i];
843
844 // Use a maximum distance for each point given by an estimation
845 // based on partition bounding boxes. The distance will be lesser
846 // than or equal to the point maximum distance.
847 double pointMaxDistance = globalMaxDistances[i];
848 for (int rank = 0; rank < m_nProcessors; ++rank) {
849 pointMaxDistance = std::min(getPartitionBox(rank).evalPointMaxDistance(point, std::numeric_limits<double>::max()), pointMaxDistance);
850 }
851
852 // Get cell distance information
853 SkdGlobalCellDistance &globalCellDistance = globalCellDistances[i];
854 int &cellRank = globalCellDistance.getRank();
855 long &cellId = globalCellDistance.getId();
856 double &cellDistance = globalCellDistance.getDistance();
857
858 // Evaluate local distance from the point
859 bool interiorCellsOnly = true;
860 nDistanceEvaluations += findPointClosestCell(point, pointMaxDistance, interiorCellsOnly, &cellId, &cellDistance);
861
862 // Set cell rank
863 if (cellId != Cell::NULL_ID) {
864 cellRank = patch.getCellOwner(cellId);
865 }
866 }
867
868 // Exchange distance information
869 MPI_Datatype globalCellDistanceDatatype = SkdGlobalCellDistance::getMPIDatatype();
870 MPI_Op globalCellDistanceMinOp = SkdGlobalCellDistance::getMPIMinOperation();
871 for (int rank = 0; rank < m_nProcessors; ++rank) {
872 SkdGlobalCellDistance *globalCellDistance = globalCellDistances.data() + globalPointsOffsets[rank];
873 if (m_rank == rank) {
874 MPI_Reduce(MPI_IN_PLACE, globalCellDistance, pointsCount[rank], globalCellDistanceDatatype, globalCellDistanceMinOp, rank, communicator);
875 } else {
876 MPI_Reduce(globalCellDistance, globalCellDistance, pointsCount[rank], globalCellDistanceDatatype, globalCellDistanceMinOp, rank, communicator);
877 }
878 }
879
880 // Update output arguments
881 for (int i = 0; i < nPoints; ++i) {
882 int globalIndex = i + globalPointsOffsets[m_rank];
883 SkdGlobalCellDistance &globalCellDistance = globalCellDistances[globalIndex];
884 globalCellDistance.exportData(ranks + i, ids + i, distances + i);
885 }
886
887 return nDistanceEvaluations;
888}
889#endif
890
891}
The PatchKernel class provides an interface for defining patches.
double getTol() const
int getCellOwner(long id) const
PatchSkdTree is the class that implements a spatial kd-tree (skd-tree) a bitpit patch.
const SkdBox & getPartitionBox(int rank) const
bool areLookupsThreadSafe() const
void clear(bool release=false)
const MPI_Comm & getCommunicator() const
const PatchKernel & getPatch() const
The SkdPatchInfo class defines a node of the skd-tree.
void updatePointClosestCell(const std::array< double, 3 > &point, bool interiorCellsOnly, long *closestId, double *closestDistance) const
void updatePointClosestCells(const std::array< double, 3 > &point, bool interiorCellsOnly, std::vector< long > *closestIds, double *closestDistance) const
static BITPIT_PUBLIC_API const std::size_t NULL_ID
The SurfaceKernel class provides an interface for defining surface patches.
long findPointClosestGlobalCell(int nPoints, const std::array< double, 3 > *points, long *ids, int *ranks, double *distances) const
void evalPointGlobalDistance(int nPoints, const std::array< double, 3 > *points, double *distances) const
void clear(bool release)
long findPointClosestCells(const std::array< double, 3 > &point, std::vector< long > &ids, double *distance) const
SurfaceSkdTree(const SurfaceKernel *patch, bool interorOnly=false)
double evalPointDistance(const std::array< double, 3 > &point) const
long findPointClosestCell(const std::array< double, 3 > &point, long *id, double *distance) const
GlobalConfigParser & root
The SkdGlobalCellDistance class allows to define a distance between a point and a cell.
static MPI_Datatype getMPIDatatype()
void exportData(int *rank, long *id, double *distance) const
--- layout: doxygen_footer ---