PeriDyno 1.0.0
Loading...
Searching...
No Matches
EdgeSet.cu
Go to the documentation of this file.
1#include "EdgeSet.h"
2#include <vector>
3#include "Array/ArrayList.h"
4
5#include <thrust/sort.h>
6
7
8namespace dyno
9{
10 template<typename TDataType>
11 EdgeSet<TDataType>::EdgeSet()
12 {
13 }
14
15 template<typename TDataType>
16 EdgeSet<TDataType>::~EdgeSet()
17 {
18 mEdges.clear();
19 mVer2Edge.clear();
20 }
21
22 __global__ void K_CountNumber(
23 DArray<uint> num,
24 DArray<TopologyModule::Edge> edges)
25 {
26 int eId = threadIdx.x + (blockIdx.x * blockDim.x);
27 if (eId >= edges.size()) return;
28
29 TopologyModule::Edge edge = edges[eId];
30
31 atomicAdd(&(num[edge[0]]), 1);
32 atomicAdd(&(num[edge[1]]), 1);
33 }
34
35 __global__ void K_StoreIds(
36 DArrayList<int> ids,
37 DArray<TopologyModule::Edge> edges)
38 {
39 int eId = threadIdx.x + (blockIdx.x * blockDim.x);
40 if (eId >= edges.size()) return;
41
42 TopologyModule::Edge edge = edges[eId];
43 int v0 = edge[0];
44 int v1 = edge[1];
45
46 ids[v0].atomicInsert(v1);
47 ids[v1].atomicInsert(v0);
48 }
49
50 template<typename TDataType>
51 void EdgeSet<TDataType>::requestPointNeighbors(DArrayList<int>& lists)
52 {
53 if (this->mCoords.isEmpty())
54 return;
55
56 DArray<uint> counts;
57 counts.resize(this->mCoords.size());
58 counts.reset();
59
60 cuExecute(mEdges.size(),
61 K_CountNumber,
62 counts,
63 mEdges);
64
65 lists.resize(counts);
66
67 cuExecute(mEdges.size(),
68 K_StoreIds,
69 lists,
70 mEdges);
71
72 counts.clear();
73 }
74
75 template<typename TDataType>
76 void EdgeSet<TDataType>::copyFrom(EdgeSet<TDataType>& edgeSet)
77 {
78 mEdges.resize(edgeSet.mEdges.size());
79 mEdges.assign(edgeSet.mEdges);
80
81 mVer2Edge.assign(edgeSet.mVer2Edge);
82
83 PointSet<TDataType>::copyFrom(edgeSet);
84 }
85
86 template<typename TDataType>
87 void EdgeSet<TDataType>::setEdges(std::vector<Edge>& edges)
88 {
89 mEdges.assign(edges);
90
91 this->tagAsChanged();
92 }
93
94 template<typename TDataType>
95 void EdgeSet<TDataType>::setEdges(DArray<Edge>& edges)
96 {
97 mEdges.resize(edges.size());
98 mEdges.assign(edges);
99
100 this->tagAsChanged();
101 }
102
103 template<typename Edge>
104 __global__ void ES_CountEdges(
105 DArray<uint> counter,
106 DArray<Edge> edges)
107 {
108 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
109 if (tId >= edges.size()) return;
110
111 Edge t = edges[tId];
112
113 atomicAdd(&counter[t[0]], 1);
114 atomicAdd(&counter[t[1]], 1);
115 }
116
117 template<typename Edge>
118 __global__ void ES_SetupEdgeIds(
119 DArrayList<int> edgeIds,
120 DArray<Edge> edges)
121 {
122 int tId = threadIdx.x + (blockIdx.x * blockDim.x);
123 if (tId >= edges.size()) return;
124
125 Edge t = edges[tId];
126
127 edgeIds[t[0]].atomicInsert(tId);
128 edgeIds[t[1]].atomicInsert(tId);
129 }
130
131 template<typename TDataType>
132 void EdgeSet<TDataType>::updateTopology()
133 {
134 this->updateEdges();
135
136 //Update the vertex to edge mapping
137 DArray<uint> counter;
138 counter.resize(this->mCoords.size());
139 counter.reset();
140
141 cuExecute(mEdges.size(),
142 ES_CountEdges,
143 counter,
144 mEdges);
145
146 mVer2Edge.resize(counter);
147
148 counter.reset();
149 cuExecute(mEdges.size(),
150 ES_SetupEdgeIds,
151 mVer2Edge,
152 mEdges);
153
154 counter.clear();
155
156 PointSet<TDataType>::updateTopology();
157 }
158
159 template<typename TDataType>
160 bool EdgeSet<TDataType>::isEmpty()
161 {
162 return mEdges.size() == 0 && PointSet<TDataType>::isEmpty();
163 }
164
165 template<typename TDataType>
166 void EdgeSet<TDataType>::clear()
167 {
168 mEdges.clear();
169 mVer2Edge.clear();
170
171 PointSet<TDataType>::clear();
172 }
173
174 DEFINE_CLASS(EdgeSet);
175}