-
Notifications
You must be signed in to change notification settings - Fork 31
Expand file tree
/
Copy pathPointSet.h
More file actions
208 lines (170 loc) · 5.26 KB
/
PointSet.h
File metadata and controls
208 lines (170 loc) · 5.26 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
#pragma once
// This is a public header. Avoid references to cuda or other external references.
#include <vector>
#include <iostream>
#include <algorithm>
#include <memory>
#include "cuda_helper.h"
#include "Common.h"
namespace cuNSearch
{
class NeighborhoodSearch;
class PointSetImplementation;
class cuNSearchDeviceData;
template<typename T, typename... Args>
std::unique_ptr<T> make_unique(Args&&... args)
{
return std::unique_ptr<T>(new T(std::forward<Args>(args)...));
}
/**
* @class PointSet.
* Represents a set of points in three-dimensional space.
*/
class PointSet
{
struct NeighborSet
{
//Pinned memory
uint NeighborCountAllocationSize;
uint ParticleCountAllocationSize;
uint *Counts;
uint *Offsets;
uint *Neighbors;
//#ifdef GPU_NEIGHBORHOOD_SEARCH
// Same data as in pinned memory. Used to avoid unnecessary data copies in the SPliSHSPlaSH-GPU implementation
uint *d_Counts;
uint *d_Offsets;
uint *d_Neighbors;
//#endif
NeighborSet()
{
NeighborCountAllocationSize = 0u;
ParticleCountAllocationSize = 0u;
Counts = nullptr;
Offsets = nullptr;
Neighbors = nullptr;
d_Counts = nullptr;
d_Offsets = nullptr;
d_Neighbors = nullptr;
}
~NeighborSet()
{
CudaHelper::CudaFree(d_Counts);
CudaHelper::CudaFree(d_Offsets);
CudaHelper::CudaFree(d_Neighbors);
}
};
public:
///**
//* Copy constructor.
//*/
PointSet(PointSet const& other);
~PointSet();
//Define descructor in cpp file to allow unique_ptr to incomplete type.
//https://stackoverflow.com/questions/9954518/stdunique-ptr-with-an-incomplete-type-wont-compile
/**
* Returns the number of neighbors of point i in the given point set.
* @param i Point index.
* @returns Number of points neighboring point i in point set point_set.
*/
inline std::size_t n_neighbors(unsigned int point_set, unsigned int i) const
{
return neighbors[point_set].Counts[i];
}
/**
* Fetches id pair of kth neighbor of point i in the given point set.
* @param point_set Point set index of other point set where neighbors have been searched.
* @param i Point index for which the neighbor id should be returned.
* @param k Represents kth neighbor of point i.
* @returns Index of neighboring point i in point set point_set.
*/
inline unsigned int neighbor(unsigned int point_set, unsigned int i, unsigned int k) const
{
//Return index of the k-th neighbor to point i (of the given point set)
const auto &neighborSet = neighbors[point_set];
return neighborSet.Neighbors[neighborSet.Offsets[i] + k];
}
inline uint n_neighborsets()
{
return neighbors.size();
}
inline uint* neighbor_indices(const uint i)
{
return neighbors[i].d_Neighbors;
}
inline uint* neighbor_counts(const uint i)
{
return neighbors[i].d_Counts;
}
inline uint* neighbor_offsets(const uint i)
{
return neighbors[i].d_Offsets;
}
PointSetImplementation *getPointSetImplementation()
{
return impl.get();
}
/**
* Fetches pointer to neighbors of point i in the given point set.
* @param point_set Point set index of other point set where neighbors have been searched.
* @param i Point index for which the neighbor id should be returned.
* @returns Pointer to ids of neighboring points of i in point set point_set.
*/
inline unsigned int * neighbor_list(unsigned int point_set, unsigned int i) const
{
//Return index of the k-th neighbor to point i (of the given point set)
const auto &neighborSet = neighbors[point_set];
return &neighborSet.Neighbors[neighborSet.Offsets[i]];
}
/**
* @returns the number of points contained in the point set.
*/
std::size_t n_points() const { return m_n; }
/*
* Returns true, if the point locations may be updated by the user.
**/
bool is_dynamic() const { return m_dynamic; }
/**
* If true is passed, the point positions may be altered by the user.
*/
void set_dynamic(bool v) { m_dynamic = v; }
Real const* GetPoints() { return m_x; }
/**
* Return the user data which can be attached to a point set.
*/
void *get_user_data() { return m_user_data; }
/**
* Reorders an array according to a previously generated sort table by invocation of the method
* "z_sort" of class "NeighborhoodSearch". Please note that the method "z_sort" of class
* "Neighborhood search" has to be called beforehand.
*/
template <typename T>
void sort_field(T* lst) const;
private:
friend NeighborhoodSearch;
friend cuNSearchDeviceData;
// Implementation and cuda data are hidden in the PointSetImplementation class to avoid unnecessary dependencies in public headers.
std::unique_ptr<PointSetImplementation> impl;
PointSet(Real const* x, std::size_t n, bool dynamic, void *user_data = nullptr);
void resize(Real const* x, std::size_t n);
Real const* point(unsigned int i) const { return &m_x[3*i]; }
Real const* m_x; //positions of the points
std::size_t m_n; //# of points in the set
bool m_dynamic; //if false the points do not move and the hash values do not change
void *m_user_data;
std::vector<uint> sortIndices;
std::vector<NeighborSet> neighbors;
};
template <typename T>
void PointSet::sort_field(T* lst) const
{
std::vector<T> tmp(lst, lst + sortIndices.size());
std::transform(sortIndices.begin(), sortIndices.end(),
//#ifdef _MSC_VER
// stdext::unchecked_array_iterator<T*>(lst),
//#else
lst,
//#endif
[&](int i) { return tmp[i]; });
}
}