PeriDyno 1.0.0
Loading...
Searching...
No Matches
SparseOctree.h
Go to the documentation of this file.
1#pragma once
2#include "DataTypes.h"
3#include "TopologyConstants.h"
5#include "Array/Array.h"
6
7namespace dyno {
8
9 typedef unsigned short OcIndex;
10 typedef unsigned long long int OcKey;
11 typedef unsigned short Level;
12
13 typedef typename ::dyno::TAlignedBox3D<Real> AABB;
14
15#define MAX_LEVEL 15
16#define DIM 3
17
19 DYN_FUNC void RecoverFromMortonCode(OcKey key, Level& l, OcIndex& x, OcIndex& y, OcIndex& z);
20
22 {
23 public:
24 DYN_FUNC OctreeNode();
25 DYN_FUNC OctreeNode(OcKey key);
26
36 DYN_FUNC OctreeNode(Level l, OcIndex x, OcIndex y, OcIndex z);
37
38 DYN_FUNC bool operator== (const OcKey k) const;
39 DYN_FUNC bool operator!= (const OcKey k) const;
40
46 DYN_FUNC bool operator>= (const OctreeNode&) const;
47
48 DYN_FUNC bool operator> (const OctreeNode&) const;
49 DYN_FUNC bool operator<= (const OctreeNode&) const;
50 DYN_FUNC bool operator< (const OctreeNode&) const;
51 DYN_FUNC bool operator== (const OctreeNode&) const;
52
53 // DYN_FUNC MortonCode3D& operator= (const MortonCode3D &);
54 // DYN_FUNC MortonCode3D operator= (const MortonCode3D &) const;
55
56 DYN_FUNC inline bool isContainedIn(const OctreeNode&) const;
57 DYN_FUNC inline bool isContainedStrictlyIn(const OctreeNode&) const;
58
59 DYN_FUNC void getCoord(Level& l, OcIndex& x, OcIndex& y, OcIndex& z);
60
61
63
64 DYN_FUNC inline OcKey key() const { return m_key; }
65 DYN_FUNC inline Level level() const { return m_level; }
66
67 DYN_FUNC inline void setDataIndex(int id) { m_data_loc = id; }
68 DYN_FUNC inline int getDataIndex() { return m_data_loc; }
69
70 DYN_FUNC inline void setStartIndex(int id) { m_start_loc = id; }
71 DYN_FUNC inline int getStartIndex() { return m_start_loc; }
72
73 DYN_FUNC inline void setFirstChildIndex(int id) { m_first_child_loc = id; }
74 DYN_FUNC inline int getFirstChildIndex() { return m_first_child_loc; }
75
76 DYN_FUNC inline void setDataSize(int n) { m_data_size = n; }
77 DYN_FUNC inline int getDataSize() { return m_data_size; }
78
79 DYN_FUNC inline bool isValid() { return m_key > 0; }
80
81 DYN_FUNC inline bool isEmpty() { return m_data_size == 0; }
82
83 public:
86
88
90 int m_data_size = 0;
91
93
95
96 bool m_bCopy = false;
97
98 int childs[8];
99 };
100
101 struct NodeCmp
102 {
103 DYN_FUNC bool operator()(const OctreeNode& A, const OctreeNode& B)
104 {
105 return A > B;
106 }
107 };
108
109 template<typename TDataType>
111 {
112 public:
113 typedef typename TDataType::Real Real;
114 typedef typename TDataType::Coord Coord;
115
118
123 void release();
124
125 void setSpace(Coord lo, Real h, Real L);
126
127 void construct(const DArray<Coord>& points, Real radius);
128 void construct(const DArray<AABB>& aabb);
129
130 void construct(const DArray<OctreeNode>& nodes);
131
132 int getLevelMax() { return m_level_max; }
133
134
136
137 GPU_FUNC Level requestLevelNumber(const AABB box);
138
139 GPU_FUNC int requestIntersectionNumber(const AABB box);
140 GPU_FUNC void reqeustIntersectionIds(int* ids, const AABB box);
141
142 GPU_FUNC int requestIntersectionNumberFromLevel(const AABB box, int level);
143 GPU_FUNC int requestIntersectionNumberFromLevel(const AABB box, AABB* data, int level);
144 GPU_FUNC void reqeustIntersectionIdsFromLevel(int* ids, const AABB box, int level);
145 GPU_FUNC void reqeustIntersectionIdsFromLevel(int* ids, const AABB box, AABB* data, int level);
146
148 GPU_FUNC void reqeustIntersectionIdsFromBottom(int* ids, const AABB box);
149
150 GPU_FUNC int requestIntersectionNumberFromBottom(const AABB box, AABB* data);
151 GPU_FUNC void reqeustIntersectionIdsFromBottom(int* ids, const AABB box, AABB* data);
152
153 public:
154
157
158 private:
159 GPU_FUNC int requestIntersectionNumber(const OcKey key, const Level l);
160 GPU_FUNC int requestIntersectionNumber(const OcKey key, const Level l, const AABB box, AABB* data);
161 GPU_FUNC void reqeustIntersectionIds(int* ids, int& shift, const OcKey key, const Level l);
162 GPU_FUNC void reqeustIntersectionIds(int* ids, int& shift, const OcKey key, const Level l, const AABB box, AABB* data);
163
164
165 private:
171
174
176
179
186 };
187}
DYN_FUNC OctreeNode(OcKey key)
DYN_FUNC OcKey key() const
DYN_FUNC bool operator==(const OcKey k) const
DYN_FUNC OctreeNode leastCommonAncestor(const OctreeNode &) const
DYN_FUNC OctreeNode(Level l, OcIndex x, OcIndex y, OcIndex z)
DYN_FUNC int getFirstChildIndex()
DYN_FUNC Level level() const
DYN_FUNC void setDataSize(int n)
DYN_FUNC bool isContainedIn(const OctreeNode &) const
DYN_FUNC bool isContainedStrictlyIn(const OctreeNode &) const
DYN_FUNC bool operator>(const OctreeNode &) const
DYN_FUNC int getStartIndex()
DYN_FUNC bool operator>=(const OctreeNode &) const
Octree will be traversed in the post-ordered fashion.
DYN_FUNC int getDataIndex()
DYN_FUNC bool operator!=(const OcKey k) const
DYN_FUNC void setDataIndex(int id)
DYN_FUNC bool isEmpty()
DYN_FUNC bool operator<(const OctreeNode &) const
DYN_FUNC void setFirstChildIndex(int id)
DYN_FUNC int getDataSize()
DYN_FUNC void setStartIndex(int id)
DYN_FUNC bool isValid()
DYN_FUNC void getCoord(Level &l, OcIndex &x, OcIndex &y, OcIndex &z)
DYN_FUNC OctreeNode()
DYN_FUNC bool operator<=(const OctreeNode &) const
GPU_FUNC int requestIntersectionNumberFromLevel(const AABB box, int level)
DArray< int > duplicates_count
GPU_FUNC void reqeustIntersectionIdsFromBottom(int *ids, const AABB box, AABB *data)
void construct(const DArray< AABB > &aabb)
DArray< OctreeNode > m_post_ordered_nodes
GPU_FUNC int requestIntersectionNumberFromBottom(const AABB box, AABB *data)
void release()
Call release() to release allocated memory explicitly, do not call this function from the decontructo...
void printPostOrderedTree()
CPU_FUNC OctreeNode queryNode(Level l, OcIndex x, OcIndex y, OcIndex z)
void construct(const DArray< OctreeNode > &nodes)
GPU_FUNC void reqeustIntersectionIds(int *ids, int &shift, const OcKey key, const Level l)
GPU_FUNC Level requestLevelNumber(const AABB box)
TDataType::Coord Coord
DArray< OctreeNode > nonRepeatNodes_cpy
DArray< OctreeNode > m_all_nodes
GPU_FUNC void reqeustIntersectionIds(int *ids, const AABB box)
GPU_FUNC int requestIntersectionNumber(const OcKey key, const Level l)
int m_level_max
levels are numbered from 0 to m_level_max;
GPU_FUNC void reqeustIntersectionIds(int *ids, int &shift, const OcKey key, const Level l, const AABB box, AABB *data)
GPU_FUNC int requestIntersectionNumber(const OcKey key, const Level l, const AABB box, AABB *data)
GPU_FUNC int requestIntersectionNumber(const AABB box)
void setSpace(Coord lo, Real h, Real L)
GPU_FUNC int requestIntersectionNumberFromLevel(const AABB box, AABB *data, int level)
GPU_FUNC void reqeustIntersectionIdsFromLevel(int *ids, const AABB box, int level)
void construct(const DArray< Coord > &points, Real radius)
GPU_FUNC void reqeustIntersectionIdsFromBottom(int *ids, const AABB box)
DArray< OctreeNode > aux_nodes
DArray< int > node_count
DArray< OctreeNode > node_buffer
TDataType::Real Real
GPU_FUNC void reqeustIntersectionIdsFromLevel(int *ids, const AABB box, AABB *data, int level)
GPU_FUNC int requestIntersectionNumberFromBottom(const AABB box)
DArray< int > data_count
This is an implementation of AdditiveCCD based on peridyno.
Definition Array.h:25
Array< T, DeviceType::GPU > DArray
Definition Array.inl:89
static DYN_FUNC void RecoverFromMortonCode(OcKey key, OcIndex &x, OcIndex &y, OcIndex &z)
Definition VoxelOctree.h:43
constexpr int EMPTY
::dyno::TAlignedBox3D< Real > AABB
static DYN_FUNC OcKey CalculateMortonCode(OcIndex x, OcIndex y, OcIndex z)
Definition VoxelOctree.h:32
unsigned long long int OcKey
DYN_FUNC bool operator()(const OctreeNode &A, const OctreeNode &B)