树结构(基于根据 Burtscher 和 Pingali 的研究)
本文将介绍广泛用于 SPH 代码和 N 体模拟中的树结构(Tree Structure)。这种数据结构主要应用于以下两个核心问题:
1. 粒子领域搜索(Neighbor Search)
在 SPH(光滑粒子流体力学)模拟中,每个粒子需要在核尺度 \( h \) 范围内查找邻居粒子,以便计算密度、压强梯度、粘性等物理量。树结构能够加速邻域搜索,尤其适用于粒子分布高度非均匀的情形。
2. 自引力计算(Self-Gravity Computation)
在引力主导的粒子系统(如星系模拟、星体碰撞)中,粒子间存在万有引力作用。直接计算所有粒子对的引力开销为 \( \mathcal{O}(N^2) \),不可接受。基于树的近似方法(如 Barnes-Hut 算法)可将计算复杂度降至 \( \mathcal{O}(N \log N) \),同时保持较高精度。
接下来的章节将分别介绍树结构在上述两个问题中的构建方法、搜索策略和性能优化。
关于传统的 链表法(Linked-List) 粒子领域搜索,请参考我另一篇博文: 👉 使用链表进行 SPH 邻域搜索
实现步骤
根据相关文献,每次执行自引力计算或粒子邻域搜索时,通常需要以下四个步骤:
确定粒子空间范围
统计所有粒子的空间边界,获取 \( x_{\min}, x_{\max}, y_{\min}, y_{\max}, z_{\min}, z_{\max} \),用于初始化树结构的根节点或空间划分范围。
构建树结构
将粒子递归划分到空间树节点中,常用的数据结构包括八叉树(Octree)或 KD 树。每个叶子节点包含若干粒子或达到最小划分条件。
粒子空间排序
对粒子进行 Morton 编码(Z-order curve)或 Hilbert 曲线编码,并按照空间位置排序,便于缓存一致性和后续并行处理。
树遍历计算
遍历树结构:
- 若执行 自引力计算,使用 Barnes-Hut 近似规则判断是否聚合节点质量;
- 若执行 邻域搜索,在每个节点中判断与查询粒子的距离是否小于核尺度 \( h \),从而筛选可能邻居。
确定粒子空间范围
直接使用规约算法求最大最小值
构建树结构
由于 GPU 上无法高效实现指针式链表结构,我们使用数组来模拟树的链接关系(例如子节点指针)。假设有一个整型数组 child,其长度远大于粒子总数 numParticles,并初始化为 -1,表示所有节点尚未使用。
区分叶子节点与根节点
在树结构中,如何区分整型数组 child 的某一个位置存储的是 叶子节点(粒子) 还是 根节点 是一个关键问题。
标识方法:
-1:表示该位置为空,未被占用。-2:表示该位置已被锁定,当前线程正在使用该位置(通常用于进行原子操作)。-3:表示该位置是一个节点,每个节点一定会有子树(可能是子节点,也可能是粒子,粒子的子节点一定为空:-1)。0 到 numParticles - 1:表示 粒子,每个位置对应一个粒子。
显然,如果第 i 个位置被锁定,那么:
- 在 三维 情况下,
8*i+1+0 到 7(即一个八岔树及其子树等)都被锁定; - 在 二维 情况下,
4*i+1+0 到 3(即一个四岔树及其子树等)都被锁定; - 在 一维 情况下,
i 及其相邻部分也会被锁定。
注意,任意索引i的子节点计算方法为:8*i+1+0,与传统的八叉树计算方法不同,这是因为我的代码中,0节点保存的是最大的根节点信息。
– 这意味着,如果一个位置被锁定,其他线程将无法访问该位置的子树或其子树的子树,确保了并行计算中节点及其相关子节点的安全。
使用原子操作同步线程
当多个线程并发地尝试访问同一个子节点槽位时,可以使用 CUDA 提供的原子操作 atomicCAS(Compare And Swap)来进行线程间同步。
以下是原子操作的代码示例:
1
| int re = atomicCAS(&child[p], -1, -2); // 尝试占用 child[p]
|
该操作的含义是:
- 若
child[p] == -1:表示该槽位尚未被占用,当前线程成功将其原子地设置为 -2,表示“锁定中”或“准备插入”; - 若
child[p] != -1:说明该槽位已被其他线程占用或已被插入,当前线程需退出或重试; - 返回值
re:表示操作前的旧值。如果 re == -1,则说明当前线程成功锁定了该节点。如果 re != -1 ,则 atomicCAS 函数发现比较不成立,直接返回了旧值,也不会替换-2而打乱child数组内容。
写到这里,树结构数组child的构建就很容易了。下面是伪代码
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
| #define childListIndex(nodeIdx, childNum) ((nodeIdx) * TREETYPE + 1 + (childNum))
#define EMPTY -1
#define LOCKED -2
#define TRUE 1
#define FALSE 0
#define NODE -3
//RealType4 是 float4 或 double4 的别名,我这样写是为了区分单双精度计算,众所周知,消费级显卡的双精度计算能力比较差。。
//粒子的w分量存放了质量,节点的w分量存放了节点半径(正方体的边长的一半)。
//目前还没有添加节点质心的计算逻辑,后续会加,用于计算粒子的引力,参考:
//A hierarchical O(N log N) force-calculation algorithm 本文发表在nature上,顶礼膜拜
__global__ void buildTreeKernel(SPHState *deviceP, treeData *tree)
{
// 为避免混淆,明确区分粒子和节点的位置数组
RealType4 *particlePositions = deviceP->positions;
volatile RealType4 *nodePositions = tree->positions;
RealType4 *nodeRoot = tree->nodeRoot;
volatile int *childList = tree->childList;
int childListLength = tree->childListLength;
int numParticles = tree->numParticles;
int inc = blockDim.x * gridDim.x;
int i = threadIdx.x + blockIdx.x * blockDim.x;
int k;
int childIndex, child;
int lockedIndex;
RealType x, y, z;
RealType rootRadius = nodeRoot->w;
RealType rootX = nodeRoot->x;
RealType rootY = nodeRoot->y;
RealType rootZ = nodeRoot->z;
// 添加跟踪当前节点位置的变量
RealType currentX, currentY, currentZ, currentR;
int depth = 0;
int isNewParticle = TRUE;
int currentNodeIndex;
bool isInsert;
while (i < numParticles)
{
isInsert = false;
while (!isInsert)
{
RealType4 pos = particlePositions[i];
depth = 0;
x = pos.x;
y = pos.y;
z = pos.z;
// 开始于根节点(索引0)
currentNodeIndex = 0;
currentX = rootX;
currentY = rootY;
currentZ = rootZ;
currentR = rootRadius;
childIndex = 0;
if (x > currentX)
childIndex = 1;
if (y > currentY)
childIndex += 2;
if (z > currentZ)
childIndex += 4;
// 跟随路径到叶节点
currentNodeIndex = childListIndex(currentNodeIndex, childIndex);
currentR *= 0.5;
currentX += ((childIndex & 1) ? currentR : -currentR);
currentY += ((childIndex & 2) ? currentR : -currentR);
currentZ += ((childIndex & 4) ? currentR : -currentR);
child = childList[currentNodeIndex];
depth++;
//下面这个while循环是为了寻找叶子节点
while (child == NODE)
{
// 确定在新节点中的子节点索引
childIndex = 0;
if (x > currentX)
childIndex = 1;
if (y > currentY)
childIndex += 2;
if (z > currentZ)
childIndex += 4;
// 跟随路径到叶节点
currentNodeIndex = childListIndex(currentNodeIndex, childIndex);
currentR *= 0.5;
currentX += ((childIndex & 1) ? currentR : -currentR);
currentY += ((childIndex & 2) ? currentR : -currentR);
currentZ += ((childIndex & 4) ? currentR : -currentR);
child = childList[currentNodeIndex];
depth++;
}
// 插入粒子到当前节点的子节点
//三种情况:
//1.当前叶子节点被占用,那么重试
//2.当前叶子节点为空,这是最简单的情况,直接插入即可
//3.当前叶子节点被粒子占用,本线程读取old节点信息,对他们两个节点进行细分,直到他们被分属到不同的象限
if (child != LOCKED)
{
lockedIndex = currentNodeIndex;
if (child == atomicCAS((int *)&childList[lockedIndex], child, LOCKED))
{
if (child == EMPTY)
{
// 直接插入粒子
childList[lockedIndex] = i;
isInsert = true;
}
else
{
//此处处理两个节点细分的情形
。。。。。。
//这个循环尝试细分,直到他们俩被分开到不同的象限
do
{
// 确定已存在粒子在新节点中的位置
int childNewIndex = 0;
if (oldParPos.x > currentX)
childNewIndex = 1;
if (oldParPos.y > currentY)
childNewIndex += 2;
if (oldParPos.z > currentZ)
childNewIndex += 4;
// 确定当前粒子在新节点中的位置
int currentNewIndex = 0;
if (x > currentX)
currentNewIndex = 1;
if (y > currentY)
currentNewIndex += 2;
if (z > currentZ)
currentNewIndex += 4;
if (childNewIndex != currentNewIndex)
{
// 两个粒子在不同子节点,可以插入
childList[childListIndex(currentNodeIndex, childNewIndex)] = child;
childList[childListIndex(currentNodeIndex, currentNewIndex)] = i;
isInsert = true;
break; // 退出循环
}
else
{
// 仍在同一子节点,需要继续细分
// 注意每次细分都会创建新节点,需要保存节点的包围盒信息 用于领域搜索
// // 写入新节点的中心和半径
nodePositions[currentNodeIndex].x = currentX;
nodePositions[currentNodeIndex].y = currentY;
nodePositions[currentNodeIndex].z = currentZ;
nodePositions[currentNodeIndex].w = currentR;
childList[currentNodeIndex] = NODE;
// 内存同步,保证所有线程都可以看到有新的节点被写入了
__threadfence();
}
} while (true);
// 确保所有子树的写入完成
__threadfence();
//释放锁
childList[lockedIndex] = NODE;
}
i += inc;
}
}
}
}
}
|
当前实现的局限性
虽然稀疏树在数据结构上直观且灵活,但在并行构建和 GPU 加速场景下,它仍然存在一些明显的局限:
- 节点访问不连续:稀疏树的节点在内存中通常分散存放,导致 GPU 在并行访问时频繁出现非连续内存读取,影响带宽利用率。
- 显存需求高:由于稀疏树存储方式不连续,GPU 在构建和查询过程中需要额外的显存来管理指针和节点结构,这使得在大规模数据下显存压力大,很容易超过显存限制。
为了解决这些问题,我们提出了 稠密存储的并行树构建方案。该方案将树节点连续存储在内存中,并结合优化的并行算法,使 GPU 能够高效地访问数据,从而显著提升构建速度和查询性能。同时,稠密存储方案能够更合理地利用显存,降低显存占用,提高处理大规模数据的能力。
更多关于实现细节和性能优化的方法,可以参考:基于稠密存储的并行树构建。
[1] Burtscher, Martin, and Keshav Pingali. “An efficient CUDA implementation of the tree-based barnes hut n-body algorithm.” GPU computing Gems Emerald edition. Morgan Kaufmann, 2011. 75-92.
title: 基于稠密存储的并行树构建实现方案
description:
date: 2025-08-15
slug: dense-tree-build
categories:
- 并行计算
- 树结构
- CUDA
draft: false
基于稠密存储的并行树构建实现方案
引言
在空间划分与邻域搜索等算法中,树形数据结构(如八叉树、四叉树、Barnes-Hut 树)是高效的加速手段。
构建这类树结构时,存储方式是影响性能与内存效率的关键因素之一。
常见的存储方式有两种:
- 稀疏存储(Sparse Storage):为每个节点预留较大索引空间,按需填充,易实现,但会浪费内存。
- 稠密存储(Dense Storage):节点存储在连续数组中,按照构建顺序紧密排列,节省内存,但需要额外的管理逻辑。
你之前的实现采用了稀疏存储,在高粒子数时内存占用明显增加。
本文将介绍如何基于稠密存储实现高效的树构建,并给出 CUDA 并行版本的实现思路。
稠密存储的基本思想
稠密存储的目标是:
- 节点数组紧凑存放,不留大块未使用空间
- 节点索引直接映射到数组下标
- 在插入新节点时,通过一个全局 maxNodeIndex 递减分配新位置
这种方法类似“倒着分配”节点空间:
1
| [粒子0] [粒子1] ... [粒子N-1] [内部节点M] [内部节点M-1] ...
|
其中:
- 0 ~ numParticles-1 区间存放叶子节点(粒子)
- numParticles ~ maxNodeIndex-1 区间存放内部节点
数据结构设计
在稠密存储中,核心数据结构通常包括:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
| // 粒子数据
struct Particles {
double *x, *y, *z; // 坐标
double *m; // 质量
double *ax, *ay, *az; // 加速度
int *depth; // 树深度
};
// 树节点信息
struct Tree {
volatile double *x, *y, *z; // 节点中心
volatile double *m; // 节点半径(或质量等)
int *childList; // 存储子节点索引
int maxNodeIndex; // 当前可分配的最大索引(递减)
};
|
关键点:childList 是一维数组,通过 childListIndex(nodeIndex, childSlot) 映射到节点的第几个子节点位置。这样存储方式天然紧凑。
构建流程详解
稠密存储的构建逻辑如下:
初始化根节点
- 根节点索引为
numNodes-1(数组末尾) - 保存中心坐标和半径
遍历粒子(并行)
- 每个线程处理多个粒子,步长为
blockDim.x * gridDim.x - 缓存粒子坐标,加速比较
下行查找插入位置
- 从根节点开始,判断粒子属于哪个子象限(八叉树中 0~7)
- 如果子节点是内部节点,继续下行
- 如果子节点是空的,直接插入粒子
- 如果子节点是另一个粒子,创建新的内部节点
创建新内部节点
- 使用
atomicSub 从 maxNodeIndex 分配新的节点索引 - 计算新节点的中心和半径
- 将已有粒子与新粒子分别插入到不同的子槽中
- 如果两者仍落在同一个槽内,继续细分
解锁与同步
- 通过
atomicCAS 实现插入位置的原子锁定 - 使用
__threadfence() 确保内存可见性
核心代码片段
下面是简化版的稠密存储节点分配逻辑:
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
| if (child == atomicCAS(&childList[lockedIndex], child, LOCKED)) {
if (child == EMPTY) {
childList[lockedIndex] = particleIndex;
} else {
// 分配新节点
newNodeIndex = atomicSub(&maxNodeIndex, 1) - 1;
// 设置新节点中心与半径
px[newNodeIndex] = currentX - r + dx;
py[newNodeIndex] = currentY - r + dy;
pz[newNodeIndex] = currentZ - r + dz;
pm[newNodeIndex] = r * 0.5;
// 初始化子节点
for (int k = 0; k < numChildren; k++) {
childList[childListIndex(newNodeIndex, k)] = EMPTY;
}
// 将旧粒子和新粒子分别放入不同槽
childList[childListIndex(newNodeIndex, oldChildSlot)] = oldParticle;
childList[childListIndex(newNodeIndex, newChildSlot)] = newParticle;
__threadfence();
childList[lockedIndex] = newNodeIndex;
}
}
|
并发与同步
由于多线程同时构建树,必须保证以下两点:
原子操作
atomicCAS(Compare And Swap)防止多个线程同时插入同一位置atomicSub 分配新节点索引
内存同步
__threadfence() 确保其他线程能看到已更新的节点数据
注意:稠密存储的节点数组在分配时是倒着使用的,所以不会和粒子索引冲突。
性能分析与优化建议
优点:
- 内存占用大幅减少(只存实际存在的节点)
- 索引紧凑,缓存命中率高
- 遍历效率提升
缺点:
- 节点分配依赖
atomicSub,在极高并发下可能成为瓶颈 - 实现复杂度高于稀疏存储
优化方向:
- 批量分配节点索引,减少
atomicSub 次数 - 合并锁与节点写入步骤,减少原子操作冲突
- 在共享内存中缓存部分子节点
总结与应用场景
稠密存储特别适合:
- 粒子数量大、空间分布均匀的模拟(如 SPH、N-body)
- 对内存占用敏感的 GPU 应用
- 需要频繁重建树结构的实时计算
与稀疏存储相比,稠密存储在 GPU 环境下通常能获得更高的性能与更低的内存消耗,尤其是在节点数量接近粒子数量时优势明显。
参考实现:本文的构建思路源自 CUDA 并行 Barnes-Hut 树构建的常见模式,并结合了你的原始代码进行稠密化处理。