xref: /libCEED/backends/hip-shared/ceed-hip-shared-basis.c (revision dd64fc8452c2d35c954858232143719e6bb2e61d)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <stdbool.h>
12 #include <stddef.h>
13 #include <hip/hip_runtime.h>
14 
15 #include "../hip/ceed-hip-common.h"
16 #include "../hip/ceed-hip-compile.h"
17 #include "ceed-hip-shared.h"
18 
19 //------------------------------------------------------------------------------
20 // Compute a block size based on required minimum threads
21 //------------------------------------------------------------------------------
22 static CeedInt ComputeBlockSizeFromRequirement(const CeedInt required) {
23   CeedInt maxSize     = 1024;  // Max total threads per block
24   CeedInt currentSize = 64;    // Start with one group
25 
26   while (currentSize < maxSize) {
27     if (currentSize > required) break;
28     else currentSize = currentSize * 2;
29   }
30   return currentSize;
31 }
32 
33 //------------------------------------------------------------------------------
34 // Compute required thread block sizes for basis kernels given P, Q, dim, and
35 // num_comp (num_comp not currently used, but may be again in other basis
36 // parallelization options)
37 //------------------------------------------------------------------------------
38 static int ComputeBasisThreadBlockSizes(const CeedInt dim, const CeedInt P_1d, const CeedInt Q_1d, const CeedInt num_comp, CeedInt *block_sizes) {
39   // Note that this will use the same block sizes for all dimensions when compiling,
40   // but as each basis object is defined for a particular dimension, we will never
41   // call any kernels except the ones for the dimension for which we have computed the
42   // block sizes.
43   const CeedInt thread_1d = CeedIntMax(P_1d, Q_1d);
44 
45   switch (dim) {
46     case 1: {
47       // Interp kernels:
48       block_sizes[0] = 256;
49 
50       // Grad kernels:
51       block_sizes[1] = 256;
52 
53       // Weight kernels:
54       block_sizes[2] = 256;
55     } break;
56     case 2: {
57       // Interp kernels:
58       CeedInt required = thread_1d * thread_1d;
59 
60       block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
61 
62       // Grad kernels: currently use same required minimum threads
63       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
64 
65       // Weight kernels:
66       required       = CeedIntMax(64, Q_1d * Q_1d);
67       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
68 
69     } break;
70     case 3: {
71       // Interp kernels:
72       CeedInt required = thread_1d * thread_1d;
73 
74       block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
75 
76       // Grad kernels: currently use same required minimum threads
77       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
78 
79       // Weight kernels:
80       required       = Q_1d * Q_1d * Q_1d;
81       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
82     }
83   }
84 
85   return CEED_ERROR_SUCCESS;
86 }
87 
88 //------------------------------------------------------------------------------
89 // Apply basis
90 //------------------------------------------------------------------------------
91 int CeedBasisApplyTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
92                                     CeedVector v) {
93   Ceed                  ceed;
94   Ceed_Hip             *ceed_Hip;
95   CeedInt               dim, num_comp;
96   const CeedScalar     *d_u;
97   CeedScalar           *d_v;
98   CeedBasis_Hip_shared *data;
99 
100   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
101   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
102   CeedCallBackend(CeedBasisGetData(basis, &data));
103   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
104   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
105 
106   // Read vectors
107   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
108   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
109   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
110 
111   // Apply basis operation
112   switch (eval_mode) {
113     case CEED_EVAL_INTERP: {
114       CeedInt P_1d, Q_1d;
115       CeedInt block_size = data->block_sizes[0];
116 
117       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
118       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
119       CeedInt thread_1d     = CeedIntMax(Q_1d, P_1d);
120       void   *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
121 
122       if (dim == 1) {
123         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
124         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
125         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
126         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
127 
128         if (t_mode == CEED_TRANSPOSE) {
129           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
130         } else {
131           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
132         }
133       } else if (dim == 2) {
134         // Check if required threads is small enough to do multiple elems
135         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
136         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
137         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
138 
139         if (t_mode == CEED_TRANSPOSE) {
140           CeedCallBackend(
141               CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
142         } else {
143           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
144         }
145       } else if (dim == 3) {
146         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
147         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
148         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
149 
150         if (t_mode == CEED_TRANSPOSE) {
151           CeedCallBackend(
152               CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
153         } else {
154           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
155         }
156       }
157     } break;
158     case CEED_EVAL_GRAD: {
159       CeedInt P_1d, Q_1d;
160       CeedInt block_size = data->block_sizes[1];
161 
162       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
163       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
164       CeedInt     thread_1d = CeedIntMax(Q_1d, P_1d);
165       CeedScalar *d_grad_1d = data->d_grad_1d;
166 
167       if (data->d_collo_grad_1d) {
168         d_grad_1d = data->d_collo_grad_1d;
169       }
170       void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v};
171       if (dim == 1) {
172         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
173         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
174         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
175         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
176 
177         if (t_mode == CEED_TRANSPOSE) {
178           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
179         } else {
180           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
181         }
182       } else if (dim == 2) {
183         // Check if required threads is small enough to do multiple elems
184         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
185         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
186         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
187 
188         if (t_mode == CEED_TRANSPOSE) {
189           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
190         } else {
191           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
192         }
193       } else if (dim == 3) {
194         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
195         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
196         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
197 
198         if (t_mode == CEED_TRANSPOSE) {
199           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
200         } else {
201           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
202         }
203       }
204     } break;
205     case CEED_EVAL_WEIGHT: {
206       CeedInt Q_1d;
207       CeedInt block_size = data->block_sizes[2];
208 
209       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
210       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
211 
212       if (dim == 1) {
213         const CeedInt opt_elems       = block_size / Q_1d;
214         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
215         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
216 
217         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
218       } else if (dim == 2) {
219         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
220         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
221         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
222 
223         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
224       } else if (dim == 3) {
225         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
226         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
227         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
228 
229         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
230       }
231     } break;
232     // LCOV_EXCL_START
233     // Evaluate the divergence to/from the quadrature points
234     case CEED_EVAL_DIV:
235       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
236     // Evaluate the curl to/from the quadrature points
237     case CEED_EVAL_CURL:
238       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
239     // Take no action, BasisApply should not have been called
240     case CEED_EVAL_NONE:
241       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
242       // LCOV_EXCL_STOP
243   }
244 
245   // Restore vectors
246   if (eval_mode != CEED_EVAL_WEIGHT) {
247     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
248   }
249   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
250   return CEED_ERROR_SUCCESS;
251 }
252 
253 //------------------------------------------------------------------------------
254 // Destroy basis
255 //------------------------------------------------------------------------------
256 static int CeedBasisDestroy_Hip_shared(CeedBasis basis) {
257   Ceed                  ceed;
258   CeedBasis_Hip_shared *data;
259 
260   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
261   CeedCallBackend(CeedBasisGetData(basis, &data));
262   CeedCallHip(ceed, hipModuleUnload(data->module));
263   CeedCallHip(ceed, hipFree(data->d_q_weight_1d));
264   CeedCallHip(ceed, hipFree(data->d_interp_1d));
265   CeedCallHip(ceed, hipFree(data->d_grad_1d));
266   CeedCallHip(ceed, hipFree(data->d_collo_grad_1d));
267   CeedCallBackend(CeedFree(&data));
268   return CEED_ERROR_SUCCESS;
269 }
270 
271 //------------------------------------------------------------------------------
272 // Create tensor basis
273 //------------------------------------------------------------------------------
274 int CeedBasisCreateTensorH1_Hip_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
275                                        const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
276   Ceed                  ceed;
277   char                 *basis_kernel_path, *basis_kernel_source;
278   CeedInt               num_comp;
279   const CeedInt         q_bytes = Q_1d * sizeof(CeedScalar);
280   const CeedInt         i_bytes = q_bytes * P_1d;
281   CeedBasis_Hip_shared *data;
282 
283   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
284   CeedCallBackend(CeedCalloc(1, &data));
285 
286   // Copy basis data to GPU
287   CeedCallHip(ceed, hipMalloc((void **)&data->d_q_weight_1d, q_bytes));
288   CeedCallHip(ceed, hipMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, hipMemcpyHostToDevice));
289   CeedCallHip(ceed, hipMalloc((void **)&data->d_interp_1d, i_bytes));
290   CeedCallHip(ceed, hipMemcpy(data->d_interp_1d, interp_1d, i_bytes, hipMemcpyHostToDevice));
291   CeedCallHip(ceed, hipMalloc((void **)&data->d_grad_1d, i_bytes));
292   CeedCallHip(ceed, hipMemcpy(data->d_grad_1d, grad_1d, i_bytes, hipMemcpyHostToDevice));
293 
294   // Compute collocated gradient and copy to GPU
295   data->d_collo_grad_1d    = NULL;
296   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
297 
298   if (has_collocated_grad) {
299     CeedScalar *collo_grad_1d;
300 
301     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
302     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
303     CeedCallHip(ceed, hipMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
304     CeedCallHip(ceed, hipMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, hipMemcpyHostToDevice));
305     CeedCallBackend(CeedFree(&collo_grad_1d));
306   }
307 
308   // Set number of threads per block for basis kernels
309   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
310   CeedCallBackend(ComputeBasisThreadBlockSizes(dim, P_1d, Q_1d, num_comp, data->block_sizes));
311 
312   // Compile basis kernels
313   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor.h", &basis_kernel_path));
314   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
315   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
316   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
317   CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->module, 11, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
318                                   CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
319                                   "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_INTERP_BLOCK_SIZE", data->block_sizes[0], "BASIS_GRAD_BLOCK_SIZE",
320                                   data->block_sizes[1], "BASIS_WEIGHT_BLOCK_SIZE", data->block_sizes[2], "BASIS_HAS_COLLOCATED_GRAD",
321                                   has_collocated_grad));
322   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Interp", &data->Interp));
323   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
324   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Grad", &data->Grad));
325   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTranspose", &data->GradTranspose));
326   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Weight", &data->Weight));
327   CeedCallBackend(CeedFree(&basis_kernel_path));
328   CeedCallBackend(CeedFree(&basis_kernel_source));
329 
330   CeedCallBackend(CeedBasisSetData(basis, data));
331 
332   // Register backend functions
333   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Hip_shared));
334   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Hip_shared));
335   return CEED_ERROR_SUCCESS;
336 }
337 
338 //------------------------------------------------------------------------------
339