forked from NVIDIA/cuda-samples
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmain.cu
325 lines (267 loc) · 11.8 KB
/
main.cu
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
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
/* Copyright (c) 2021, NVIDIA CORPORATION. All rights reserved.
*
* Redistribution and use in source and binary forms, with or without
* modification, are permitted provided that the following conditions
* are met:
* * Redistributions of source code must retain the above copyright
* notice, this list of conditions and the following disclaimer.
* * Redistributions in binary form must reproduce the above copyright
* notice, this list of conditions and the following disclaimer in the
* documentation and/or other materials provided with the distribution.
* * Neither the name of NVIDIA CORPORATION nor the names of its
* contributors may be used to endorse or promote products derived
* from this software without specific prior written permission.
*
* THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS ``AS IS'' AND ANY
* EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
* IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
* PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT OWNER OR
* CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
* EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
* PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
* PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY
* OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
* (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
* OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*/
/* Computation of eigenvalues of symmetric, tridiagonal matrix using
* bisection.
*/
// includes, system
#include <stdlib.h>
#include <stdio.h>
#include <string.h>
#include <math.h>
#include <float.h>
#include <assert.h>
// includes, project
#include <helper_functions.h>
#include <helper_cuda.h>
#include "config.h"
#include "structs.h"
#include "matlab.h"
#include "util.h"
#include "gerschgorin.h"
#include "bisect_small.cuh"
#include "bisect_large.cuh"
////////////////////////////////////////////////////////////////////////////////
// declaration, forward
bool runTest(int argc, char **argv);
////////////////////////////////////////////////////////////////////////////////
// Program main
////////////////////////////////////////////////////////////////////////////////
int main(int argc, char **argv) {
bool bQAResults = false;
printf("Starting eigenvalues\n");
bQAResults = runTest(argc, argv);
printf("Test %s\n", bQAResults ? "Succeeded!" : "Failed!");
exit(bQAResults ? EXIT_SUCCESS : EXIT_FAILURE);
}
////////////////////////////////////////////////////////////////////////////////
//! Initialize the input data to the algorithm
//! @param input handles to the input data
//! @param exec_path path where executable is run (argv[0])
//! @param mat_size size of the matrix
//! @param user_defined 1 if the matrix size has been requested by the user,
//! 0 if the default size
////////////////////////////////////////////////////////////////////////////////
void initInputData(InputData &input, char *exec_path,
const unsigned int mat_size,
const unsigned int user_defined) {
// allocate memory
input.a = (float *)malloc(sizeof(float) * mat_size);
input.b = (float *)malloc(sizeof(float) * mat_size);
if (1 == user_defined) {
// initialize diagonal and superdiagonal entries with random values
srand(278217421);
// srand( clock());
for (unsigned int i = 0; i < mat_size; ++i) {
input.a[i] = (float)(2.0 * (((double)rand() / (double)RAND_MAX) - 0.5));
input.b[i] = (float)(2.0 * (((double)rand() / (double)RAND_MAX) - 0.5));
}
// the first element of s is used as padding on the device (thus the
// whole vector is copied to the device but the kernels are launched
// with (s+1) as start address
input.b[0] = 0.0f;
} else {
// read default matrix
unsigned int input_data_size = mat_size;
char *diag_path = sdkFindFilePath("diagonal.dat", exec_path);
assert(NULL != diag_path);
sdkReadFile(diag_path, &(input.a), &input_data_size, false);
char *sdiag_path = sdkFindFilePath("superdiagonal.dat", exec_path);
assert(NULL != sdiag_path);
sdkReadFile(sdiag_path, &(input.b), &input_data_size, false);
free(diag_path);
free(sdiag_path);
}
// allocate device memory for input
checkCudaErrors(cudaMalloc((void **)&(input.g_a), sizeof(float) * mat_size));
checkCudaErrors(
cudaMalloc((void **)&(input.g_b_raw), sizeof(float) * mat_size));
// copy data to device
checkCudaErrors(cudaMemcpy(input.g_a, input.a, sizeof(float) * mat_size,
cudaMemcpyHostToDevice));
checkCudaErrors(cudaMemcpy(input.g_b_raw, input.b, sizeof(float) * mat_size,
cudaMemcpyHostToDevice));
input.g_b = input.g_b_raw + 1;
}
////////////////////////////////////////////////////////////////////////////////
//! Clean up input data, in particular allocated memory
//! @param input handles to the input data
////////////////////////////////////////////////////////////////////////////////
void cleanupInputData(InputData &input) {
freePtr(input.a);
freePtr(input.b);
checkCudaErrors(cudaFree(input.g_a));
input.g_a = NULL;
checkCudaErrors(cudaFree(input.g_b_raw));
input.g_b_raw = NULL;
input.g_b = NULL;
}
////////////////////////////////////////////////////////////////////////////////
//! Check if a specific matrix size has to be used
//! @param argc number of command line arguments (from main(argc, argv)
//! @param argv pointers to command line arguments (from main(argc, argv)
//! @param matrix_size size of matrix, updated if specific size specified on
//! command line
////////////////////////////////////////////////////////////////////////////////
void getMatrixSize(int argc, char **argv, unsigned int &mat_size,
unsigned int &user_defined) {
int temp = -1;
if (checkCmdLineFlag(argc, (const char **)argv, "matrix-size")) {
temp = getCmdLineArgumentInt(argc, (const char **)argv, "matrix-size");
}
if (temp > 0) {
mat_size = (unsigned int)temp;
// data type short is used in the kernel
assert(mat_size < (1 << 16));
// mat_size should be large than 2
assert(mat_size >= 2);
user_defined = 1;
}
printf("Matrix size: %i x %i\n", mat_size, mat_size);
}
////////////////////////////////////////////////////////////////////////////////
//! Check if a specific precision of the eigenvalue has to be obtained
//! @param argc number of command line arguments (from main(argc, argv)
//! @param argv pointers to command line arguments (from main(argc, argv)
//! @param iters_timing numbers of iterations for timing, updated if a
//! specific number is specified on the command line
//! @param user_defined 1 if the precision has been requested by the user,
//! 0 if the default size
////////////////////////////////////////////////////////////////////////////////
void getPrecision(int argc, char **argv, float &precision,
unsigned int &user_defined) {
float temp = -1.0f;
if (checkCmdLineFlag(argc, (const char **)argv, "precision")) {
temp = getCmdLineArgumentFloat(argc, (const char **)argv, "precision");
printf("Precision is between [0.001, 0.000001]\n");
}
if (temp > 1e-6 && temp <= 0.001) {
precision = temp;
user_defined = 1;
}
printf("Precision: %f\n", precision);
}
////////////////////////////////////////////////////////////////////////////////
//! Check if a particular number of iterations for timings has to be used
//! @param argc number of command line arguments (from main(argc, argv)
//! @param argv pointers to command line arguments (from main(argc, argv)
//! @param iters_timing number of timing iterations, updated if user
//! specific value
////////////////////////////////////////////////////////////////////////////////
void getItersTiming(int argc, char **argv, unsigned int &iters_timing) {
int temp = -1;
if (checkCmdLineFlag(argc, (const char **)argv, "iters-timing")) {
temp = getCmdLineArgumentInt(argc, (const char **)argv, "iters-timing");
}
if (temp > 0) {
iters_timing = temp;
}
printf("Iterations to be timed: %i\n", iters_timing);
}
////////////////////////////////////////////////////////////////////////////////
//! Check if a particular filename has to be used for the file where the result
//! is stored
//! @param argc number of command line arguments (from main(argc, argv)
//! @param argv pointers to command line arguments (from main(argc, argv)
//! @param filename filename of result file, updated if user specified
//! filename
////////////////////////////////////////////////////////////////////////////////
void getResultFilename(int argc, char **argv, char *&filename) {
char *temp = NULL;
getCmdLineArgumentString(argc, (const char **)argv, "filename-result", &temp);
if (NULL != temp) {
filename = (char *)malloc(sizeof(char) * strlen(temp));
strcpy(filename, temp);
free(temp);
}
printf("Result filename: '%s'\n", filename);
}
////////////////////////////////////////////////////////////////////////////////
//! Run a simple test for CUDA
////////////////////////////////////////////////////////////////////////////////
bool runTest(int argc, char **argv) {
bool bCompareResult = false;
findCudaDevice(argc, (const char **)argv);
StopWatchInterface *timer = NULL;
StopWatchInterface *timer_total = NULL;
sdkCreateTimer(&timer);
sdkCreateTimer(&timer_total);
// default
unsigned int mat_size = 2048;
// flag if the matrix size is due to explicit user request
unsigned int user_defined = 0;
// desired precision of eigenvalues
float precision = 0.00001f;
unsigned int iters_timing = 100;
char *result_file = (char *)"eigenvalues.dat";
// check if there is a command line request for the matrix size
getMatrixSize(argc, argv, mat_size, user_defined);
// check if user requested specific precision
getPrecision(argc, argv, precision, user_defined);
// check if user requested specific number of iterations for timing
getItersTiming(argc, argv, iters_timing);
// file name for result file
getResultFilename(argc, argv, result_file);
// set up input
InputData input;
initInputData(input, argv[0], mat_size, user_defined);
// compute Gerschgorin interval
float lg = FLT_MAX;
float ug = -FLT_MAX;
computeGerschgorin(input.a, input.b + 1, mat_size, lg, ug);
printf("Gerschgorin interval: %f / %f\n", lg, ug);
// two kernels, for small matrices a lot of overhead can be avoided
if (mat_size <= MAX_SMALL_MATRIX) {
// initialize memory for result
ResultDataSmall result;
initResultSmallMatrix(result, mat_size);
// run the kernel
computeEigenvaluesSmallMatrix(input, result, mat_size, lg, ug, precision,
iters_timing);
// get the result from the device and do some sanity checks,
// save the result
processResultSmallMatrix(input, result, mat_size, result_file);
// clean up
cleanupResultSmallMatrix(result);
printf("User requests non-default argument(s), skipping self-check!\n");
bCompareResult = true;
} else {
// initialize memory for result
ResultDataLarge result;
initResultDataLargeMatrix(result, mat_size);
// run the kernel
computeEigenvaluesLargeMatrix(input, result, mat_size, precision, lg, ug,
iters_timing);
// get the result from the device and do some sanity checks
// save the result if user specified matrix size
bCompareResult = processResultDataLargeMatrix(
input, result, mat_size, result_file, user_defined, argv[0]);
// cleanup
cleanupResultDataLargeMatrix(result);
}
cleanupInputData(input);
return bCompareResult;
}