BlosSOM
Interactive dimensionality reduction on large datasets (EmbedSOM and FLOWER combined)
embedsom_cuda.cpp
Go to the documentation of this file.
1/*
2The MIT License
3
4Copyright (c) 2021 Martin Krulis
5 Mirek Kratochvil
6 Sona Molnarova
7
8Permission is hereby granted, free of charge,
9to any person obtaining a copy of this software and
10associated documentation files (the "Software"), to
11deal in the Software without restriction, including
12without limitation the rights to use, copy, modify,
13merge, publish, distribute, sublicense, and/or sell
14copies of the Software, and to permit persons to whom
15the Software is furnished to do so,
16subject to the following conditions:
17
18The above copyright notice and this permission notice
19shall be included in all copies or substantial portions of the Software.
20
21THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
22EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES
23OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
24IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR
25ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
26TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
27SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
28*/
29
30#include "embedsom_cuda.h"
31
32#include "cuda_runtime.h"
33
34void
36 size_t g,
37 size_t d,
38 float boost,
39 size_t k,
40 float adjust,
41 const float *hidim_points,
42 const float *hidim_landmarks,
43 const float *lodim_landmarks,
44 float *lodim_points)
45{
46
47 /* first make sure all sizes are acceptable and there's actual sensible
48 * work that can be done */
49 if (!n)
50 return;
51 if (d < 2 || d > 256)
52 return; // TODO try d=1; d>256 should warn
53 if (k > g)
54 k = g;
55 if (k < 3)
56 return;
57
58 /* make sure the buffers are sufficiently big. The allocations are kept as
59 * such: |data| = d*n |points| = 2*n |lm_hi| = d*g |lm_lo| = 2*g |topk| =
60 * adjusted_k()*n
61 *
62 * Because `n` and `g` is highly variable in this usecase, we _only_
63 * _increase_ the buffer sizes. Reallocation in case the sizes grow totally
64 * out of limits would help too (feel free to add into the conditions).
65 */
66
67 CUCH(cudaSetDevice(0));
68
69#define size_check(tgt, cur, buffer, type) \
70 if (tgt > cur) { \
71 if (buffer) \
72 cudaFree(buffer); \
73 CUCH(cudaMalloc(&buffer, tgt * sizeof(type))); \
74 cur = tgt; \
75 }
76
77 size_check(d * n, ndata, data, float);
78 size_check(2 * n, npoints, points, float);
79 size_check(d * g, nlm_hi, lm_hi, float);
80 size_check(2 * g, nlm_lo, lm_lo, float);
81
82 size_t adjusted_k = std::min(g, k + 1);
83 size_check(adjusted_k * n, nknns, knns, knn_entry<float>);
84 CUCH(cudaGetLastError());
85#undef size_check
86
87 /* all is safely allocated, transfer the data! */
88 CUCH(cudaMemcpy(
89 data, hidim_points, d * n * sizeof(float), cudaMemcpyHostToDevice));
90 CUCH(cudaMemcpy(
91 lm_hi, hidim_landmarks, d * g * sizeof(float), cudaMemcpyHostToDevice));
92 CUCH(cudaMemcpy(
93 lm_lo, lodim_landmarks, 2 * g * sizeof(float), cudaMemcpyHostToDevice));
94 CUCH(cudaGetLastError());
95
96 /* run the main embedding */
97 runKNNKernel(d, n, g, adjusted_k);
98 runProjectionKernel(d, n, g, k, boost, adjust);
99
100 /* get the results */
101 CUCH(cudaMemcpy(
102 lodim_points, points, 2 * n * sizeof(float), cudaMemcpyDeviceToHost));
103 CUCH(cudaGetLastError());
104}
105
107{
108 if (data)
109 cudaFree(data);
110 if (points)
111 cudaFree(points);
112 if (lm_hi)
113 cudaFree(lm_hi);
114 if (lm_lo)
115 cudaFree(lm_lo);
116 if (knns)
117 cudaFree(knns);
118}
#define size_check(tgt, cur, buffer, type)
#define CUCH(status)
Macro wrapper for CUDA calls checking.
void runKNNKernel(size_t d, size_t n, size_t g, size_t adjusted_k)
Execute the CUDA kNN within the context.
void runProjectionKernel(size_t d, size_t n, size_t g, size_t k, float boost, float adjust)
Execute the CUDA projection kernel atop the context.
void run(size_t n, size_t g, size_t d, float boost, size_t k, float adjust, const float *hidim_points, const float *hidim_landmarks, const float *lodim_landmarks, float *lodim_points)
Run EmbedSOM on the given data.
knn_entry< float > * knns
Definition: embedsom_cuda.h:54