source: CIVL/examples/cuda/matMult1.cu@ af389e5

1.23 2.0 main test-branch
Last change on this file since af389e5 was 92d17822, checked in by Andre Marianiello <andre.marianiello@…>, 11 years ago

Converted matMult1.cu example to work with symbolic inputs, and added a test for it to Cuda2CIVLTransformTest.java

git-svn-id: svn://vsl.cis.udel.edu/civl/trunk@1966 fb995dde-84ed-4084-dfe6-e5aef3e2452c

  • Property mode set to 100644
File size: 8.9 KB
Line 
1/***********************************************************************
2* FILENAME: MM.cu
3* Matrix Multiplication
4* Matrix operands have row-major order.
5*
6* C = A * B
7* Multiplies two square matrices (NxN * NxN).
8* Matrix values have type double.
9*
10* A simple CUDA program has a basic workflow:
11* 1) Initialize matrix operands as double-precision arrays on host (CPU).
12* 2) Copy operands from host memory to GPU memory.
13* 3) Apply matrix operaton to operands on GPU
14* 4) Copy result from GPU memory to host memory.
15*
16*
17* CUDA C Programming Guide Version 4.2 (3.2.3, p.22):
18* http://developer.download.nvidia.com/compute/DevZone/docs/html/C/doc/CUDA_C_Programming_Guide.pdf
19*
20* MM with linearized matrix operands:
21* http://www.hpcwire.com/hpcwire/2008-10-08/compilers_and_more_programming_gpus_today.html
22*
23*************************************************************************/
24// online source: https://www.rcac.purdue.edu/userinfo/resources/carter/compile/MM.cu
25
26#include <stdio.h>
27#include <stdlib.h>
28#include "cuda.h"
29
30#ifdef _CIVL
31$input int N;
32$input int TILE_WIDTH;
33#else
34#define N 1024 /* size of square matrix */
35#define TILE_WIDTH 16
36#endif
37
38
39/* MM kernel using global (not shared) memory. */
40__global__
41void myMM_global (const double * const A, const double * const B, double *C, int width) {
42
43 /* Get row and column from block and thread IDs */
44 int row = (blockDim.y*blockIdx.y) + threadIdx.y;
45 int col = (blockDim.x*blockIdx.x) + threadIdx.x;
46
47 /* Initialize result of one element which one thread computes. */
48 double result=0.0;
49
50 /* Compute one element of the matrix product. */
51 for (int i = 0; i < width; ++i)
52 result += A[row*width + i] * B[i*width + col];
53
54 /* Store the result of one matrix element in matrix C. */
55 C[row * width + col] = result;
56}
57
58
59/* MM kernel using shared memory. */
60__global__
61void myMM_shared (const double * const A, const double * const B, double* C, int width) {
62 __shared__ double A_shared[TILE_WIDTH][TILE_WIDTH];
63 __shared__ double B_shared[TILE_WIDTH][TILE_WIDTH];
64
65 int bx = blockIdx.x; int by = blockIdx.y;
66 int tx = threadIdx.x; int ty = threadIdx.y;
67
68 /* Identify the row and column of the C element to work on. */
69 int row = by * TILE_WIDTH + ty;
70 int col = bx * TILE_WIDTH + tx;
71
72 double result = 0.0;
73
74 /* Loop over the A and B tiles required to compute the C element. */
75 for (int phase = 0; phase < width/TILE_WIDTH; ++phase) {
76 /* Shared effort: loading of A and B tiles into shared memory. */
77 A_shared[ty][tx] = A[row*width + (phase*TILE_WIDTH + tx)];
78 B_shared[ty][tx] = B[col + (phase*TILE_WIDTH + ty)*width];
79 __syncthreads();
80
81 for (int k = 0; k < TILE_WIDTH; ++k)
82 result += A_shared[ty][k] * B_shared[k][tx];
83 __syncthreads();
84
85 }
86 C[row*width+col] = result;
87}
88
89
90/************************************************************************/
91/************************************************************************/
92/************************************************************************/
93
94
95int main (int argc, char** argv) {
96
97#ifdef _CIVL
98 $assume argc == 2;
99 $assume atoi(argv[1]) == 0;
100#endif
101
102 /* Set device based on input from command line */
103 if (argc > 1) {
104 if (cudaSetDevice(atoi(argv[1])) != cudaSuccess) {
105 int num_devices;
106 cudaGetDeviceCount(&num_devices);
107 fprintf(stderr, "Error initializing device %s,\
108 device value must be 0-%d\n", argv[1], (num_devices-1));
109 return 0;
110 }
111 } else {
112 fprintf(stderr, "Usage: %s gpu_device\n", argv[0]);
113 return 0;
114 }
115
116 /* Declare CPU arrays. */
117#ifdef _CIVL
118 $input double A[N*N];
119 $input double B[N*N];
120#else
121 double A[N*N],B[N*N];
122#endif
123 double C[N*N]; /* linearized CPU double arrays */
124 int r,c;
125
126 /* Declare GPU arrays. */
127 double *G_A,*G_B,*G_C; /* linearized GPU double arrays */
128 size_t size_a,size_b,size_c; /* size of linearized array in bytes */
129 size_a = size_b = size_c = N*N;
130
131 /* Setup a clock. */
132 cudaEvent_t start, stop;
133 float CPU_elapsedtime, GPU_global_elapsedtime, GPU_shared_elapsedtime;
134 cudaEventCreate(&start);
135 cudaEventCreate(&stop);
136
137
138
139
140 /* 1) Initialize matrix operands as double-precision arrays on host (CPU). */
141#ifndef _CIVL
142 for (r=0;r<N;++r)
143 for (c=0;c<N;++c) {
144 A[r*N+c] = 1.0;
145 B[r*N+c] = 1.0;
146 }
147#endif
148
149
150/*-----------------------------------------------------------------------*/
151
152 /* MM on a CPU. */
153 cudaEventRecord(start,0);
154 for (int r = 0; r < N; ++r )
155 for (int c = 0; c < N; ++c )
156 for (int k = 0; k < N; ++k )
157 C[r*N+c] += A[r*N+c] * B[k*N+c];
158 cudaEventRecord(stop,0);
159 cudaEventSynchronize(stop);
160 cudaEventElapsedTime(&CPU_elapsedtime,start,stop);
161 printf(" speedup\n");
162 printf(" -------\n");
163 printf("Elapsed time in CPU: %7.1f milliseconds\n", CPU_elapsedtime);
164/*-----------------------------------------------------------------------*/
165
166 /* MM on Global Memory of GPGPU. */
167 cudaEventRecord(start,0);
168
169 /* 2) Copy operands from CPU memory to GPGPU memory. */
170 cudaMalloc((void**)&G_A,size_a*sizeof(double)); /* alloc A in GPGPU */
171 cudaMalloc((void**)&G_B,size_b*sizeof(double)); /* alloc B in GPGPU */
172 cudaMalloc((void**)&G_C,size_c*sizeof(double)); /* alloc C in GPGPU */
173 cudaMemcpy(G_A,A,size_a*sizeof(double),cudaMemcpyHostToDevice);
174 cudaMemcpy(G_B,B,size_b*sizeof(double),cudaMemcpyHostToDevice);
175
176 /* 3) Apply matrix operation to operands on GPGPU */
177 /* There is no partial final block in this example. */
178 dim3 block = {TILE_WIDTH,TILE_WIDTH,1}; /* using a 2D block: 16,16,1 */
179 dim3 grid = {N/TILE_WIDTH,N/TILE_WIDTH,1}; /* as many 16x16-thread blocks as needed: */
180 myMM_global<<< grid,block >>>(G_A,G_B,G_C,N); /* grid(16,16,1) */
181
182 /* 4) Copy result from GPGPU memory to CPU memory. */
183 cudaMemcpy(C,G_C,size_c*sizeof(double),cudaMemcpyDeviceToHost);
184
185 /* Deallocate memory on GPGPU. */
186 cudaFree(G_A);
187 cudaFree(G_B);
188 cudaFree(G_C);
189
190 cudaEventRecord(stop,0);
191 cudaEventSynchronize(stop);
192 cudaEventElapsedTime(&GPU_global_elapsedtime,start,stop);
193 printf("Elapsed time in GPU (global memory): %7.1f milliseconds %5.1f\n",
194 GPU_global_elapsedtime,CPU_elapsedtime/GPU_global_elapsedtime);
195//*
196 printf("\nGLOBAL MEMORY:\n");
197 for (r=0;r<N;++r)
198 for (c=0;c<N;++c) {
199 printf("%2d,%2d %g\n", r,c,C[r*N+c]);
200 }
201//*/
202/*-----------------------------------------------------------------------*/
203
204 /* MM on Shared Memory of GPGPU. */
205 cudaEventRecord(start,0);
206
207 /* 2) Copy operands from CPU memory to GPGPU memory. */
208 cudaMalloc((void**)&G_A,size_a*sizeof(double)); /* alloc A in GPGPU */
209 cudaMalloc((void**)&G_B,size_b*sizeof(double)); /* alloc B in GPGPU */
210 cudaMalloc((void**)&G_C,size_c*sizeof(double)); /* alloc C in GPGPU */
211 cudaMemcpy(G_A,A,size_a*sizeof(double),cudaMemcpyHostToDevice);
212 cudaMemcpy(G_B,B,size_b*sizeof(double),cudaMemcpyHostToDevice);
213
214 /* 3) Apply matrix operation to operands on GPGPU */
215 /* There is not partial final block in this example. */
216 /* Use the same grid and block from the previous case. */
217 myMM_shared<<< grid,block >>>(G_A,G_B,G_C,N);
218
219 /* 4) Copy result from GPGPU memory to CPU memory. */
220 cudaMemcpy(C,G_C,size_c*sizeof(double),cudaMemcpyDeviceToHost);
221
222 /* Deallocate memory on GPGPU. */
223 cudaFree(G_A);
224 cudaFree(G_B);
225 cudaFree(G_C);
226
227 cudaEventRecord(stop,0);
228 cudaEventSynchronize(stop);
229 cudaEventElapsedTime(&GPU_shared_elapsedtime,start,stop);
230 printf("Elapsed time in GPU (shared memory): %7.1f milliseconds %5.1f\n",
231 GPU_shared_elapsedtime,CPU_elapsedtime/GPU_shared_elapsedtime);
232//*
233 printf("\nSHARED MEMORY:\n");
234 for (r=0;r<N;++r)
235 for (c=0;c<N;++c) {
236 printf("%2d,%2d %g\n", r,c,C[r*N+c]);
237 }
238//*/
239/*-----------------------------------------------------------------------*/
240
241 /* Deallocate the clock. */
242 cudaEventDestroy(start);
243 cudaEventDestroy(stop);
244
245 return 0;
246}
247
Note: See TracBrowser for help on using the repository browser.