PeriDyno 1.0.0
Loading...
Searching...
No Matches
SurfaceInteraction.cu
Go to the documentation of this file.
1#include "SurfaceInteraction.h"
2#include <thrust/sort.h>
3#include <iostream>
4#include <OrbitCamera.h>
5
6namespace dyno
7{
8 __global__ void SI_SurfaceInitializeArray(
9 DArray<int> intersected)
10 {
11 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
12 if (pId >= intersected.size()) return;
13
14 intersected[pId] = 0;
15 }
16
17 __global__ void SI_SurfaceMergeIntersectedIndexOR(
18 DArray<int> intersected1,
19 DArray<int> intersected2,
20 DArray<int> outIntersected,
21 DArray<int> outUnintersected)
22 {
23 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
24 if (pId >= intersected1.size()) return;
25
26 if (intersected1[pId] == 0 && intersected2[pId] == 0)
27 outIntersected[pId] = 0;
28 else
29 outIntersected[pId] = 1;
30
31 outUnintersected[pId] = outIntersected[pId] == 1 ? 0 : 1;
32 }
33
34 __global__ void SI_SurfaceMergeIntersectedIndexXOR(
35 DArray<int> intersected1,
36 DArray<int> intersected2,
37 DArray<int> outIntersected,
38 DArray<int> outUnintersected)
39 {
40 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
41 if (pId >= intersected1.size()) return;
42
43 if (intersected1[pId] == intersected2[pId])
44 outIntersected[pId] = 0;
45 else
46 outIntersected[pId] = 1;
47
48 outUnintersected[pId] = outIntersected[pId] == 1 ? 0 : 1;
49 }
50
51 __global__ void SI_SurfaceMergeIntersectedIndexC(
52 DArray<int> intersected1,
53 DArray<int> intersected2,
54 DArray<int> outIntersected,
55 DArray<int> outUnintersected)
56 {
57 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
58 if (pId >= intersected1.size()) return;
59
60 if (intersected2[pId] == 1)
61 outIntersected[pId] = 0;
62 else
63 outIntersected[pId] = intersected1[pId];
64
65 outUnintersected[pId] = outIntersected[pId] == 1 ? 0 : 1;
66 }
67
68 template<typename TDataType>
69 SurfaceInteraction<TDataType>::SurfaceInteraction()
70 {
71 this->ray1 = TRay3D<Real>();
72 this->ray2 = TRay3D<Real>();
73 this->isPressed = false;
74
75 this->outOtherTriangleSet()->setDataPtr(std::make_shared<TriangleSet<TDataType>>());
76 this->outOtherTriangleSet()->getDataPtr()->getTriangles().resize(0);
77 this->outSelectedTriangleSet()->setDataPtr(std::make_shared<TriangleSet<TDataType>>());
78 this->outSelectedTriangleSet()->getDataPtr()->getTriangles().resize(0);
79 }
80
81 template<typename TDataType>
82 void SurfaceInteraction<TDataType>::onEvent(PMouseEvent event)
83 {
84 if (!event.altKeyPressed()) {
85 if (camera == nullptr)
86 {
87 this->camera = event.camera;
88 }
89 this->varToggleMultiSelect()->setValue(false);
90 if (event.shiftKeyPressed() || event.controlKeyPressed())
91 {
92 this->varToggleMultiSelect()->setValue(true);
93 if (event.shiftKeyPressed() && !event.controlKeyPressed())
94 {
95 this->varMultiSelectionType()->getDataPtr()->setCurrentKey(0);
96 }
97 else if (!event.shiftKeyPressed() && event.controlKeyPressed())
98 {
99 this->varMultiSelectionType()->getDataPtr()->setCurrentKey(1);
100 }
101 else if (event.shiftKeyPressed() && event.controlKeyPressed())
102 {
103 this->varMultiSelectionType()->getDataPtr()->setCurrentKey(2);;
104 }
105 }
106 if (event.actionType == AT_PRESS)
107 {
108 this->camera = event.camera;
109 this->isPressed = true;
110 //printf("Mouse pressed: Origin: %f %f %f; Direction: %f %f %f \n", event.ray.origin.x, event.ray.origin.y, event.ray.origin.z, event.ray.direction.x, event.ray.direction.y, event.ray.direction.z);
111 this->ray1.origin = event.ray.origin;
112 this->ray1.direction = event.ray.direction;
113 this->x1 = event.x;
114 this->y1 = event.y;
115 if (this->varSurfacePickingType()->getValue() == PickingTypeSelection::Both || this->varSurfacePickingType()->getValue() == PickingTypeSelection::Click)
116 {
117 this->calcIntersectClick();
118 this->printInfoClick();
119 }
120 }
121 else if (event.actionType == AT_RELEASE)
122 {
123 this->isPressed = false;
124 //printf("Mouse released: Origin: %f %f %f; Direction: %f %f %f \n", event.ray.origin.x, event.ray.origin.y, event.ray.origin.z, event.ray.direction.x, event.ray.direction.y, event.ray.direction.z);
125 this->ray2.origin = event.ray.origin;
126 this->ray2.direction = event.ray.direction;
127 this->x2 = event.x;
128 this->y2 = event.y;
129 if (this->varToggleMultiSelect()->getValue() && this->varTogglePicker()->getValue())
130 {
131 this->mergeIndex();
132 this->printInfoDragRelease();
133 }
134 }
135 else
136 {
137 //printf("Mouse repeated: Origin: %f %f %f; Direction: %f %f %f \n", event.ray.origin.x, event.ray.origin.y, event.ray.origin.z, event.ray.direction.x, event.ray.direction.y, event.ray.direction.z);
138 if (this->isPressed)
139 {
140 this->ray2.origin = event.ray.origin;
141 this->ray2.direction = event.ray.direction;
142 this->x2 = event.x;
143 this->y2 = event.y;
144 if (abs(this->x2 - this->x1) <= 3 && abs(this->y2 - this->y1) <= 3)
145 {
146 if (this->varSurfacePickingType()->getValue() == PickingTypeSelection::Both || this->varSurfacePickingType()->getValue() == PickingTypeSelection::Click)
147 this->calcIntersectClick();
148 }
149 else
150 {
151 if (this->varSurfacePickingType()->getValue() == PickingTypeSelection::Both || this->varSurfacePickingType()->getValue() == PickingTypeSelection::Drag)
152 {
153 this->calcIntersectDrag();
154 this->printInfoDragging();
155 }
156 }
157 }
158 }
159 }
160 }
161
162 template <typename Triangle, typename Real, typename Coord>
163 __global__ void SI_CalIntersectedTrisRay(
164 DArray<Coord> points,
165 DArray<Triangle> triangles,
166 DArray<int> intersected,
167 DArray<int> unintersected,
168 DArray<Real> triDistance,
169 TRay3D<Real> mouseray)
170 {
171 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
172 if (pId >= triangles.size()) return;
173
174 TTriangle3D<Real> t = TTriangle3D<Real>(points[triangles[pId][0]], points[triangles[pId][1]], points[triangles[pId][2]]);
175 int temp = 0;
176
177 TPoint3D<Real> p;
178 temp = mouseray.intersect(t, p);
179
180 if (temp == 1)
181 {
182 intersected[pId] = 1;
183 triDistance[pId] = (mouseray.origin - p.origin).norm();
184 }
185 else
186 {
187 intersected[pId] = 0;
188 triDistance[pId] = 3.4E38;
189 }
190 unintersected[pId] = (intersected[pId] == 1 ? 0 : 1);
191 }
192
193 __global__ void SI_CalTrisNearest(
194 int min_index,
195 DArray<int> intersected,
196 DArray<int> unintersected)
197 {
198 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
199 if (pId >= intersected.size()) return;
200
201 if (intersected[pId] == 1)
202 {
203 if (pId != min_index)
204 {
205 intersected[pId] = 0;
206 unintersected[pId] = 1;
207 }
208 }
209 }
210
211 template <typename Triangle, typename Real, typename Coord>
212 __global__ void SI_CalIntersectedTrisBox(
213 DArray<Coord> points,
214 DArray<Triangle> triangles,
215 DArray<int> intersected,
216 DArray<int> unintersected,
217 TPlane3D<Real> plane13,
218 TPlane3D<Real> plane42,
219 TPlane3D<Real> plane14,
220 TPlane3D<Real> plane32,
221 TRay3D<Real> mouseray)
222 {
223 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
224 if (pId >= triangles.size()) return;
225
226 TTriangle3D<Real> t = TTriangle3D<Real>(points[triangles[pId][0]], points[triangles[pId][1]], points[triangles[pId][2]]);
227 TSegment3D<Real> s1 = TSegment3D<Real>(points[triangles[pId][0]], points[triangles[pId][1]]);
228 TSegment3D<Real> s2 = TSegment3D<Real>(points[triangles[pId][1]], points[triangles[pId][2]]);
229 TSegment3D<Real> s3 = TSegment3D<Real>(points[triangles[pId][2]], points[triangles[pId][0]]);
230
231 bool flag = false;
232 TPoint3D<Real> p;
233 bool temp11 = s1.intersect(plane13, p);
234 temp11 = temp11 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
235 if (temp11)
236 flag = true;
237 bool temp12 = s1.intersect(plane42, p);
238 temp12 = temp12 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
239 if (temp12)
240 flag = true;
241 bool temp13 = s1.intersect(plane14, p);
242 temp13 = temp13 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
243 if (temp13)
244 flag = true;
245 bool temp14 = s1.intersect(plane32, p);
246 temp14 = temp14 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
247 if (temp14)
248 flag = true;
249 bool temp21 = s2.intersect(plane13, p);
250 temp21 = temp21 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
251 if (temp21)
252 flag = true;
253 bool temp22 = s2.intersect(plane42, p);
254 temp22 = temp22 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
255 if (temp22)
256 flag = true;
257 bool temp23 = s2.intersect(plane14, p);
258 temp23 = temp23 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
259 if (temp23)
260 flag = true;
261 bool temp24 = s2.intersect(plane32, p);
262 temp24 = temp24 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
263 if (temp24)
264 flag = true;
265 bool temp31 = s3.intersect(plane13, p);
266 temp31 = temp31 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
267 if (temp31)
268 flag = true;
269 bool temp32 = s3.intersect(plane42, p);
270 temp32 = temp32 && (((p.origin - plane14.origin).dot(plane14.normal)) * ((p.origin - plane32.origin).dot(plane32.normal))) > 0;
271 if (temp32)
272 flag = true;
273 bool temp33 = s3.intersect(plane14, p);
274 temp33 = temp33 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
275 if (temp33)
276 flag = true;
277 bool temp34 = s3.intersect(plane32, p);
278 temp34 = temp34 && (((p.origin - plane13.origin).dot(plane13.normal)) * ((p.origin - plane42.origin).dot(plane42.normal))) > 0;
279 if (temp34)
280 flag = true;
281
282 for (int i = 0; i < 3; i++)
283 {
284 float temp1 = ((points[triangles[pId][i]] - plane13.origin).dot(plane13.normal)) * ((points[triangles[pId][i]] - plane42.origin).dot(plane42.normal));
285 float temp2 = ((points[triangles[pId][i]] - plane14.origin).dot(plane14.normal)) * ((points[triangles[pId][i]] - plane32.origin).dot(plane32.normal));
286 if (temp1 >= 0 && temp2 >= 0)
287 {
288 flag = flag || true;
289 break;
290 }
291 }
292
293 int temp = mouseray.intersect(t, p);
294 if (temp == 1)
295 {
296 flag = flag || true;
297 }
298
299 if (flag || intersected[pId] == 1)
300 intersected[pId] = 1;
301 else
302 intersected[pId] = 0;
303 unintersected[pId] = (intersected[pId] == 1 ? 0 : 1);
304 }
305
306 template <typename Triangle>
307 __global__ void SI_AssignOutTriangles(
308 DArray<Triangle> triangles,
309 DArray<Triangle> intersected_triangles,
310 DArray<Triangle> unintersected_triangles,
311 DArray<int> outTriangleIndex,
312 DArray<int> intersected,
313 DArray<int> unintersected,
314 DArray<int> intersected_o)
315 {
316 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
317 if (pId >= triangles.size()) return;
318
319 if (intersected_o[pId] == 1)
320 {
321 intersected_triangles[intersected[pId]] = triangles[pId];
322 outTriangleIndex[intersected[pId]] = pId;
323 }
324 else
325 {
326 unintersected_triangles[unintersected[pId]] = triangles[pId];
327 }
328 }
329
330 template <typename Triangle, typename Coord, typename Real>
331 __global__ void SI_NeighborTrisDiffuse(
332 DArray<Triangle> triangles,
333 DArray<Coord> points,
334 DArray<TopologyModule::Tri2Edg> tri2Edg,
335 DArray<TopologyModule::Edg2Tri> edg2Tri,
336 DArray<int> intersected,
337 DArray<int> unintersected,
338 Real diffusionAngle)
339 {
340 int i = threadIdx.x + blockIdx.x * blockDim.x;
341 if (i >= triangles.size()) return;
342
343 if (intersected[i] == 1)
344 {
345 TTriangle3D<Real> t0, t1;
346 t0 = TTriangle3D<Real>(points[triangles[i][0]], points[triangles[i][1]], points[triangles[i][2]]);
347 for (int ii = 0; ii < 3; ii++)
348 {
349 for (int iii = 0; iii < 2; iii++)
350 {
351 if (intersected[edg2Tri[tri2Edg[i][ii]][iii]] != 1)
352 {
353 t1 = TTriangle3D<Real>(points[triangles[edg2Tri[tri2Edg[i][ii]][iii]][0]], points[triangles[edg2Tri[tri2Edg[i][ii]][iii]][1]], points[triangles[edg2Tri[tri2Edg[i][ii]][iii]][2]]);
354 if (abs(t0.normal().dot(t1.normal())) >= cosf(diffusionAngle))
355 {
356 intersected[edg2Tri[tri2Edg[i][ii]][iii]] = 1;
357 unintersected[edg2Tri[tri2Edg[i][ii]][iii]] = 0;
358 }
359 }
360 }
361 }
362 }
363 }
364
365 template <typename Triangle, typename Coord, typename Real>
366 __global__ void SI_VisibleFilter(
367 DArray<Triangle> triangles,
368 DArray<Coord> points,
369 DArray<int> intersected,
370 DArray<int> unintersected,
371 TRay3D<Real> mouseray)
372 {
373 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
374 if (pId >= triangles.size()) return;
375
376 if (intersected[pId] == 1)
377 {
378 TTriangle3D<Real> t = TTriangle3D<Real>(points[triangles[pId][0]], points[triangles[pId][1]], points[triangles[pId][2]]);
379 if (mouseray.direction.dot(t.normal()) >= 0)
380 {
381 intersected[pId] = 0;
382 unintersected[pId] = 1;
383 }
384 }
385 }
386
387 __global__ void SI_Tri2Quad(
388 DArray<int> intersected,
389 DArray<int> unintersected)
390 {
391 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
392 if (pId >= intersected.size()) return;
393
394 if (intersected[pId] == 1)
395 {
396 if (pId % 2 == 1)
397 {
398 if (intersected[pId - 1] == 0)
399 {
400 intersected[pId - 1] = 1;
401 unintersected[pId - 1] = 0;
402 }
403 }
404 else
405 {
406
407 if (intersected[pId + 1] == 0)
408 {
409 intersected[pId + 1] = 1;
410 unintersected[pId + 1] = 0;
411 }
412 }
413 }
414 }
415
416 __global__ void SI_QuadOutput(
417 DArray<int> intersected_o,
418 DArray<int> intersected_q
419 )
420 {
421 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
422 if (pId >= intersected_q.size()) return;
423
424 intersected_q[pId] = 0;
425 if (intersected_o[pId*2] == 1||intersected_o[pId*2+1]==1)
426 {
427 intersected_q[pId] = 1;
428 }
429 }
430
431 __global__ void SI_QuadIndexOutput(
432 DArray<int> outTriangleIndex,
433 DArray<int> outQuadIndex
434 )
435 {
436 int pId = threadIdx.x + (blockIdx.x * blockDim.x);
437 if (pId >= outQuadIndex.size()) return;
438
439 outQuadIndex[pId] = outTriangleIndex[pId * 2] / 2;
440 }
441
442 __global__ void SI_InitialS2PSelected(
443 DArray<int> s2PSelected)
444 {
445 int pId = threadIdx.x + blockIdx.x * blockDim.x;
446 if (pId >= s2PSelected.size()) return;
447
448 s2PSelected[pId] = 0;
449 }
450
451 template <typename Triangle>
452 __global__ void SI_Surface2Point(
453 DArray<Triangle> triangles,
454 DArray<int> outTriangleIndex,
455 DArray<int> s2PSelected)
456 {
457 int pId = threadIdx.x + blockIdx.x * blockDim.x;
458 if (pId >= outTriangleIndex.size()) return;
459
460 s2PSelected[triangles[outTriangleIndex[pId]][0]] = 1;
461 s2PSelected[triangles[outTriangleIndex[pId]][1]] = 1;
462 s2PSelected[triangles[outTriangleIndex[pId]][2]] = 1;
463 }
464
465 __global__ void SI_s2PIndexOut(
466 DArray<int> s2PIndex,
467 DArray<int> s2PIndex_o,
468 DArray<int> s2PIndexOut)
469 {
470 int pId = threadIdx.x + blockIdx.x * blockDim.x;
471 if (pId >= s2PIndex.size()) return;
472
473 if (s2PIndex_o[pId] == 1)
474 {
475 s2PIndexOut[s2PIndex[pId]] = pId;
476 }
477 }
478
479
480 template<typename TDataType>
481 void SurfaceInteraction<TDataType>::calcSurfaceIntersectClick()
482 {
483 auto& initialTriangleSet = this->inInitialTriangleSet()->getData();
484 auto& points = initialTriangleSet.getPoints();
485 auto& triangles = initialTriangleSet.getTriangles();
486 DArray<int> intersected;
487 intersected.resize(triangles.size());
488 cuExecute(triangles.size(),
489 SI_SurfaceInitializeArray,
490 intersected
491 );
492 DArray<int> unintersected;
493 unintersected.resize(triangles.size());
494 this->tempNumT = triangles.size();
495
496 DArray<Real> triDistance;
497 triDistance.resize(triangles.size());
498
499 cuExecute(triangles.size(),
500 SI_CalIntersectedTrisRay,
501 points,
502 triangles,
503 intersected,
504 unintersected,
505 triDistance,
506 this->ray1
507 );
508
509 int min_index = thrust::min_element(thrust::device, triDistance.begin(), triDistance.begin() + triDistance.size()) - triDistance.begin();
510
511 cuExecute(intersected.size(),
512 SI_CalTrisNearest,
513 min_index,
514 intersected,
515 unintersected
516 );
517
518 if (this->varToggleFlood()->getValue())
519 {
520 auto& Tri2Edg = initialTriangleSet.getTriangle2Edge();
521 auto& Edge2Tri = initialTriangleSet.getEdge2Triangle();
522 Edge2Tri.resize(0);
523 initialTriangleSet.updateTriangle2Edge();
524 int intersected_size_t_o = 0;
525 int intersected_size_t = 1;
526 while (intersected_size_t > intersected_size_t_o && intersected_size_t < triangles.size())
527 {
528 intersected_size_t_o = intersected_size_t;
529
530 cuExecute(triangles.size(),
531 SI_NeighborTrisDiffuse,
532 triangles,
533 points,
534 Tri2Edg,
535 Edge2Tri,
536 intersected,
537 unintersected,
538 Real(this->varFloodAngle()->getValue()/180.0f*M_PI));
539
540 intersected_size_t = thrust::reduce(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), (int)0, thrust::plus<int>());
541 }
542 }
543
544 if (this->varToggleVisibleFilter()->getValue())
545 {
546 cuExecute(triangles.size(),
547 SI_VisibleFilter,
548 triangles,
549 points,
550 intersected,
551 unintersected,
552 this->ray1
553 );
554 }
555
556 if (this->varToggleQuad()->getValue())
557 {
558 cuExecute(triangles.size(),
559 SI_Tri2Quad,
560 intersected,
561 unintersected
562 )
563 }
564
565 this->tempTriIntersectedIndex.assign(intersected);
566
567 if (this->varToggleMultiSelect()->getData())
568 {
569 if (this->triIntersectedIndex.size() == 0)
570 {
571 this->triIntersectedIndex.resize(triangles.size());
572 cuExecute(triangles.size(),
573 SI_SurfaceInitializeArray,
574 this->triIntersectedIndex
575 )
576 }
577 DArray<int> outIntersected;
578 outIntersected.resize(intersected.size());
579 DArray<int> outUnintersected;
580 outUnintersected.resize(unintersected.size());
581 if (this->varMultiSelectionType()->getValue() == MultiSelectionType::OR)
582 {
583 cuExecute(triangles.size(),
584 SI_SurfaceMergeIntersectedIndexOR,
585 this->triIntersectedIndex,
586 intersected,
587 outIntersected,
588 outUnintersected
589 );
590 }
591 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::XOR)
592 {
593 cuExecute(triangles.size(),
594 SI_SurfaceMergeIntersectedIndexXOR,
595 this->triIntersectedIndex,
596 intersected,
597 outIntersected,
598 outUnintersected
599 );
600 }
601 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::C)
602 {
603 cuExecute(triangles.size(),
604 SI_SurfaceMergeIntersectedIndexC,
605 this->triIntersectedIndex,
606 intersected,
607 outIntersected,
608 outUnintersected
609 );
610 }
611 intersected.assign(outIntersected);
612 unintersected.assign(outUnintersected);
613 }
614 else
615 {
616 this->triIntersectedIndex.assign(intersected);
617 }
618
619 DArray<int> intersected_o;
620 intersected_o.assign(intersected);
621
622 int intersected_size = thrust::reduce(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), (int)0, thrust::plus<int>());
623 DArray<int> outTriangleIndex;
624 outTriangleIndex.resize(intersected_size);
625 thrust::exclusive_scan(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), intersected.begin());
626 DArray<Triangle> intersected_triangles;
627 intersected_triangles.resize(intersected_size);
628
629 int unintersected_size = thrust::reduce(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), (int)0, thrust::plus<int>());
630 thrust::exclusive_scan(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), unintersected.begin());
631 DArray<Triangle> unintersected_triangles;
632 unintersected_triangles.resize(unintersected_size);
633
634 cuExecute(triangles.size(),
635 SI_AssignOutTriangles,
636 triangles,
637 intersected_triangles,
638 unintersected_triangles,
639 outTriangleIndex,
640 intersected,
641 unintersected,
642 intersected_o
643 );
644
645 DArray<int> s2PSelected;
646 s2PSelected.resize(points.size());
647
648 cuExecute(points.size(),
649 SI_InitialS2PSelected,
650 s2PSelected
651 );
652
653 cuExecute(outTriangleIndex.size(),
654 SI_Surface2Point,
655 triangles,
656 outTriangleIndex,
657 s2PSelected
658 );
659 int s2PSelectedSize = thrust::reduce(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), (int)0, thrust::plus<int>());
660 DArray<int> s2PSelected_o;
661 s2PSelected_o.assign(s2PSelected);
662
663 thrust::exclusive_scan(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), s2PSelected.begin());
664
665 DArray<int> s2PSelectedIndex;
666 s2PSelectedIndex.resize(s2PSelectedSize);
667 cuExecute(s2PSelected.size(),
668 SI_s2PIndexOut,
669 s2PSelected,
670 s2PSelected_o,
671 s2PSelectedIndex
672 );
673
674 if (this->varToggleQuad()->getValue())
675 {
676 DArray<int> intersected_q;
677 intersected_q.resize(triangles.size() / 2);
678 cuExecute(triangles.size(),
679 SI_QuadOutput,
680 intersected_o,
681 intersected_q
682 );
683 intersected_o.assign(intersected_q);
684
685 DArray<int> outQuadIndex;
686 outQuadIndex.resize(outTriangleIndex.size()/2);
687 cuExecute(outQuadIndex.size(),
688 SI_QuadIndexOutput,
689 outTriangleIndex,
690 outQuadIndex
691 );
692 outTriangleIndex.assign(outQuadIndex);
693 }
694
695 this->tempNumS = intersected_size;
696 this->outSelectedTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
697 this->outSelectedTriangleSet()->getDataPtr()->setTriangles(intersected_triangles);
698 this->outOtherTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
699 this->outOtherTriangleSet()->getDataPtr()->setTriangles(unintersected_triangles);
700 if (this->varToggleIndexOutput()->getValue())
701 {
702 this->outTriangleIndex()->getDataPtr()->assign(outTriangleIndex);
703 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelectedIndex);
704 }
705 else
706 {
707 this->outTriangleIndex()->getDataPtr()->assign(intersected_o);
708 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelected_o);
709 }
710 }
711
712 template<typename TDataType>
713 void SurfaceInteraction<TDataType>::calcSurfaceIntersectDrag()
714 {
715 if (x1 == x2)
716 {
717 x2 += 1.0f;
718 }
719 if (y1 == y2)
720 {
721 y2 += 1.0f;
722 }
723 TRay3D<Real> ray1 = this->camera->castRayInWorldSpace((float)x1, (float)y1);
724 TRay3D<Real> ray2 = this->camera->castRayInWorldSpace((float)x2, (float)y2);
725 TRay3D<Real> ray3 = this->camera->castRayInWorldSpace((float)x1, (float)y2);
726 TRay3D<Real> ray4 = this->camera->castRayInWorldSpace((float)x2, (float)y1);
727
728 TPlane3D<Real> plane13 = TPlane3D<Real>(ray1.origin, ray1.direction.cross(ray3.direction));
729 TPlane3D<Real> plane42 = TPlane3D<Real>(ray2.origin, ray2.direction.cross(ray4.direction));
730 TPlane3D<Real> plane14 = TPlane3D<Real>(ray4.origin, ray1.direction.cross(ray4.direction));
731 TPlane3D<Real> plane32 = TPlane3D<Real>(ray3.origin, ray2.direction.cross(ray3.direction));
732
733 auto& initialTriangleSet = this->inInitialTriangleSet()->getData();
734 auto& points = initialTriangleSet.getPoints();
735 auto& triangles = initialTriangleSet.getTriangles();
736 DArray<int> intersected;
737 intersected.resize(triangles.size());
738 cuExecute(triangles.size(),
739 SI_SurfaceInitializeArray,
740 intersected
741 );
742 DArray<int> unintersected;
743 unintersected.resize(triangles.size());
744 this->tempNumT = triangles.size();
745
746 cuExecute(triangles.size(),
747 SI_CalIntersectedTrisBox,
748 points,
749 triangles,
750 intersected,
751 unintersected,
752 plane13,
753 plane42,
754 plane14,
755 plane32,
756 this->ray2
757 );
758
759 if (this->varToggleVisibleFilter()->getValue())
760 {
761 cuExecute(triangles.size(),
762 SI_VisibleFilter,
763 triangles,
764 points,
765 intersected,
766 unintersected,
767 this->ray1
768 );
769 }
770
771 if (this->varToggleQuad()->getValue())
772 {
773 cuExecute(triangles.size(),
774 SI_Tri2Quad,
775 intersected,
776 unintersected
777 )
778 }
779
780 this->tempTriIntersectedIndex.assign(intersected);
781
782 if (this->varToggleMultiSelect()->getData())
783 {
784 if (this->triIntersectedIndex.size() == 0)
785 {
786 this->triIntersectedIndex.resize(triangles.size());
787 cuExecute(triangles.size(),
788 SI_SurfaceInitializeArray,
789 this->triIntersectedIndex
790 )
791 }
792 DArray<int> outIntersected;
793 outIntersected.resize(intersected.size());
794 DArray<int> outUnintersected;
795 outUnintersected.resize(unintersected.size());
796 if (this->varMultiSelectionType()->getValue() == MultiSelectionType::OR)
797 {
798 cuExecute(triangles.size(),
799 SI_SurfaceMergeIntersectedIndexOR,
800 this->triIntersectedIndex,
801 intersected,
802 outIntersected,
803 outUnintersected
804 );
805 }
806 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::XOR)
807 {
808 cuExecute(triangles.size(),
809 SI_SurfaceMergeIntersectedIndexXOR,
810 this->triIntersectedIndex,
811 intersected,
812 outIntersected,
813 outUnintersected
814 );
815 }
816 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::C)
817 {
818 cuExecute(triangles.size(),
819 SI_SurfaceMergeIntersectedIndexC,
820 this->triIntersectedIndex,
821 intersected,
822 outIntersected,
823 outUnintersected
824 );
825 }
826 intersected.assign(outIntersected);
827 unintersected.assign(outUnintersected);
828 //this->triIntersectedIndex.assign(intersected);
829 }
830 else
831 {
832 this->triIntersectedIndex.assign(intersected);
833 }
834
835
836 DArray<int> intersected_o;
837 intersected_o.assign(intersected);
838
839 int intersected_size = thrust::reduce(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), (int)0, thrust::plus<int>());
840 DArray<int> outTriangleIndex;
841 outTriangleIndex.resize(intersected_size);
842 thrust::exclusive_scan(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), intersected.begin());
843 DArray<Triangle> intersected_triangles;
844 intersected_triangles.resize(intersected_size);
845
846 int unintersected_size = thrust::reduce(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), (int)0, thrust::plus<int>());
847 thrust::exclusive_scan(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), unintersected.begin());
848 DArray<Triangle> unintersected_triangles;
849 unintersected_triangles.resize(unintersected_size);
850
851 cuExecute(triangles.size(),
852 SI_AssignOutTriangles,
853 triangles,
854 intersected_triangles,
855 unintersected_triangles,
856 outTriangleIndex,
857 intersected,
858 unintersected,
859 intersected_o
860 );
861
862 DArray<int> s2PSelected;
863 s2PSelected.resize(points.size());
864
865 cuExecute(points.size(),
866 SI_InitialS2PSelected,
867 s2PSelected
868 );
869
870 cuExecute(outTriangleIndex.size(),
871 SI_Surface2Point,
872 triangles,
873 outTriangleIndex,
874 s2PSelected
875 );
876 int s2PSelectedSize = thrust::reduce(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), (int)0, thrust::plus<int>());
877 DArray<int> s2PSelected_o;
878 s2PSelected_o.assign(s2PSelected);
879
880 thrust::exclusive_scan(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), s2PSelected.begin());
881
882 DArray<int> s2PSelectedIndex;
883 s2PSelectedIndex.resize(s2PSelectedSize);
884 cuExecute(s2PSelected.size(),
885 SI_s2PIndexOut,
886 s2PSelected,
887 s2PSelected_o,
888 s2PSelectedIndex
889 );
890
891 if (this->varToggleQuad()->getValue())
892 {
893 DArray<int> intersected_q;
894 intersected_q.resize(triangles.size() / 2);
895 cuExecute(triangles.size(),
896 SI_QuadOutput,
897 intersected_o,
898 intersected_q
899 );
900 intersected_o.assign(intersected_q);
901
902 DArray<int> outQuadIndex;
903 outQuadIndex.resize(outTriangleIndex.size() / 2);
904 cuExecute(outQuadIndex.size(),
905 SI_QuadIndexOutput,
906 outTriangleIndex,
907 outQuadIndex
908 );
909 outTriangleIndex.assign(outQuadIndex);
910 }
911
912 this->tempNumS = intersected_size;
913 this->outSelectedTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
914 this->outSelectedTriangleSet()->getDataPtr()->setTriangles(intersected_triangles);
915 this->outOtherTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
916 this->outOtherTriangleSet()->getDataPtr()->setTriangles(unintersected_triangles);
917 if (this->varToggleIndexOutput()->getValue())
918 {
919 this->outTriangleIndex()->getDataPtr()->assign(outTriangleIndex);
920 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelectedIndex);
921 }
922 else
923 {
924 this->outTriangleIndex()->getDataPtr()->assign(intersected_o);
925 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelected_o);
926 }
927 }
928
929 template<typename TDataType>
930 void SurfaceInteraction<TDataType>::mergeIndex()
931 {
932 auto& initialTriangleSet = this->inInitialTriangleSet()->getData();
933 auto& points = initialTriangleSet.getPoints();
934 auto& triangles = initialTriangleSet.getTriangles();
935 DArray<int> intersected;
936 intersected.resize(triangles.size());
937 cuExecute(triangles.size(),
938 SI_SurfaceInitializeArray,
939 intersected
940 );
941 DArray<int> unintersected;
942 unintersected.resize(triangles.size());
943 this->tempNumT = triangles.size();
944
945 DArray<int> outIntersected;
946 outIntersected.resize(intersected.size());
947 DArray<int> outUnintersected;
948 outUnintersected.resize(unintersected.size());
949
950 if (this->varMultiSelectionType()->getValue() == MultiSelectionType::OR)
951 {
952 cuExecute(triangles.size(),
953 SI_SurfaceMergeIntersectedIndexOR,
954 this->triIntersectedIndex,
955 this->tempTriIntersectedIndex,
956 outIntersected,
957 outUnintersected
958 );
959 }
960 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::XOR)
961 {
962 cuExecute(triangles.size(),
963 SI_SurfaceMergeIntersectedIndexXOR,
964 this->triIntersectedIndex,
965 this->tempTriIntersectedIndex,
966 outIntersected,
967 outUnintersected
968 );
969 }
970 else if (this->varMultiSelectionType()->getValue() == MultiSelectionType::C)
971 {
972 cuExecute(triangles.size(),
973 SI_SurfaceMergeIntersectedIndexC,
974 this->triIntersectedIndex,
975 this->tempTriIntersectedIndex,
976 outIntersected,
977 outUnintersected
978 );
979 }
980
981 intersected.assign(outIntersected);
982 unintersected.assign(outUnintersected);
983 this->triIntersectedIndex.assign(intersected);
984
985 DArray<int> intersected_o;
986 intersected_o.assign(intersected);
987
988 int intersected_size = thrust::reduce(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), (int)0, thrust::plus<int>());
989 DArray<int> outTriangleIndex;
990 outTriangleIndex.resize(intersected_size);
991 thrust::exclusive_scan(thrust::device, intersected.begin(), intersected.begin() + intersected.size(), intersected.begin());
992 DArray<Triangle> intersected_triangles;
993 intersected_triangles.resize(intersected_size);
994
995 int unintersected_size = thrust::reduce(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), (int)0, thrust::plus<int>());
996 thrust::exclusive_scan(thrust::device, unintersected.begin(), unintersected.begin() + unintersected.size(), unintersected.begin());
997 DArray<Triangle> unintersected_triangles;
998 unintersected_triangles.resize(unintersected_size);
999
1000 cuExecute(triangles.size(),
1001 SI_AssignOutTriangles,
1002 triangles,
1003 intersected_triangles,
1004 unintersected_triangles,
1005 outTriangleIndex,
1006 intersected,
1007 unintersected,
1008 intersected_o
1009 );
1010
1011 DArray<int> s2PSelected;
1012 s2PSelected.resize(points.size());
1013
1014 cuExecute(points.size(),
1015 SI_InitialS2PSelected,
1016 s2PSelected
1017 );
1018
1019 cuExecute(outTriangleIndex.size(),
1020 SI_Surface2Point,
1021 triangles,
1022 outTriangleIndex,
1023 s2PSelected
1024 );
1025 int s2PSelectedSize = thrust::reduce(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), (int)0, thrust::plus<int>());
1026 DArray<int> s2PSelected_o;
1027 s2PSelected_o.assign(s2PSelected);
1028
1029 thrust::exclusive_scan(thrust::device, s2PSelected.begin(), s2PSelected.begin() + s2PSelected.size(), s2PSelected.begin());
1030
1031 DArray<int> s2PSelectedIndex;
1032 s2PSelectedIndex.resize(s2PSelectedSize);
1033 cuExecute(s2PSelected.size(),
1034 SI_s2PIndexOut,
1035 s2PSelected,
1036 s2PSelected_o,
1037 s2PSelectedIndex
1038 );
1039
1040 this->tempNumS = intersected_size;
1041 this->outSelectedTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
1042 this->outSelectedTriangleSet()->getDataPtr()->setTriangles(intersected_triangles);
1043 this->outOtherTriangleSet()->getDataPtr()->copyFrom(initialTriangleSet);
1044 this->outOtherTriangleSet()->getDataPtr()->setTriangles(unintersected_triangles);
1045 if (this->varToggleIndexOutput()->getValue())
1046 {
1047 this->outTriangleIndex()->getDataPtr()->assign(outTriangleIndex);
1048 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelectedIndex);
1049 }
1050 else
1051 {
1052 this->outTriangleIndex()->getDataPtr()->assign(intersected_o);
1053 this->outSur2PointIndex()->getDataPtr()->assign(s2PSelected_o);
1054 }
1055 }
1056
1057 template<typename TDataType>
1058 void SurfaceInteraction<TDataType>::printInfoClick()
1059 {
1060 std::cout << "----------surface picking: click----------" << std::endl;
1061 std::cout << "multiple picking: " << this->varToggleMultiSelect()->getValue() << std::endl;
1062 std::cout << "selected num/ total num:" << this->tempNumS << "/" << this->tempNumT << std::endl;
1063 }
1064
1065 template<typename TDataType>
1066 void SurfaceInteraction<TDataType>::printInfoDragging()
1067 {
1068 std::cout << "----------surface picking: dragging----------" << std::endl;
1069 std::cout << "multiple picking: " << this->varToggleMultiSelect()->getValue() << std::endl;
1070 std::cout << "selected num/ total num:" << this->tempNumS << "/" << this->tempNumT << std::endl;
1071 }
1072
1073 template<typename TDataType>
1074 void SurfaceInteraction<TDataType>::printInfoDragRelease()
1075 {
1076 std::cout << "----------surface picking: drag release----------" << std::endl;
1077 std::cout << "multiple picking: " << this->varToggleMultiSelect()->getValue() << std::endl;
1078 std::cout << "selected num/ total num:" << this->tempNumS << "/" << this->tempNumT << std::endl;
1079 }
1080
1081 template<typename TDataType>
1082 void SurfaceInteraction<TDataType>::calcIntersectClick()
1083 {
1084 if (this->varTogglePicker()->getData())
1085 calcSurfaceIntersectClick();
1086 }
1087
1088 template<typename TDataType>
1089 void SurfaceInteraction<TDataType>::calcIntersectDrag()
1090 {
1091 if (this->varTogglePicker()->getData())
1092 calcSurfaceIntersectDrag();
1093 }
1094
1095 DEFINE_CLASS(SurfaceInteraction);
1096}