xref: /libCEED/backends/hip-shared/ceed-hip-shared-basis.c (revision fa5adaf534a16ab6b728a4ab386a8a8b7691fe89)
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   switch (dim) {
45     case 1: {
46       // Interp kernels:
47       block_sizes[0] = 256;
48 
49       // Grad kernels:
50       block_sizes[1] = 256;
51 
52       // Weight kernels:
53       block_sizes[2] = 256;
54     } break;
55     case 2: {
56       // Interp kernels:
57       CeedInt required = thread_1d * thread_1d;
58       block_sizes[0]   = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
59 
60       // Grad kernels: currently use same required minimum threads
61       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
62 
63       // Weight kernels:
64       required       = CeedIntMax(64, Q_1d * Q_1d);
65       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
66 
67     } break;
68     case 3: {
69       // Interp kernels:
70       CeedInt required = thread_1d * thread_1d;
71       block_sizes[0]   = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
72 
73       // Grad kernels: currently use same required minimum threads
74       block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
75 
76       // Weight kernels:
77       required       = Q_1d * Q_1d * Q_1d;
78       block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required));
79     }
80   }
81 
82   return CEED_ERROR_SUCCESS;
83 }
84 
85 //------------------------------------------------------------------------------
86 // Apply basis
87 //------------------------------------------------------------------------------
88 int CeedBasisApplyTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
89                                     CeedVector v) {
90   Ceed ceed;
91   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
92   Ceed_Hip *ceed_Hip;
93   CeedCallBackend(CeedGetData(ceed, &ceed_Hip));
94   CeedBasis_Hip_shared *data;
95   CeedCallBackend(CeedBasisGetData(basis, &data));
96   CeedInt dim, num_comp;
97   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
98   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
99 
100   // Read vectors
101   const CeedScalar *d_u;
102   CeedScalar       *d_v;
103   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
104   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
105   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
106 
107   // Apply basis operation
108   switch (eval_mode) {
109     case CEED_EVAL_INTERP: {
110       CeedInt P_1d, Q_1d;
111       CeedInt block_size = data->block_sizes[0];
112       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
113       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
114       CeedInt thread_1d     = CeedIntMax(Q_1d, P_1d);
115       void   *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v};
116       if (dim == 1) {
117         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
118         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
119         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
120         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
121         if (t_mode == CEED_TRANSPOSE) {
122           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
123         } else {
124           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
125         }
126       } else if (dim == 2) {
127         // Check if required threads is small enough to do multiple elems
128         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
129         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
130         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
131         if (t_mode == CEED_TRANSPOSE) {
132           CeedCallBackend(
133               CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
134         } else {
135           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
136         }
137       } else if (dim == 3) {
138         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
139         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
140         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
141         if (t_mode == CEED_TRANSPOSE) {
142           CeedCallBackend(
143               CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
144         } else {
145           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
146         }
147       }
148     } break;
149     case CEED_EVAL_GRAD: {
150       CeedInt P_1d, Q_1d;
151       CeedInt block_size = data->block_sizes[1];
152       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
153       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
154       CeedInt     thread_1d = CeedIntMax(Q_1d, P_1d);
155       CeedScalar *d_grad_1d = data->d_grad_1d;
156       if (data->d_collo_grad_1d) {
157         d_grad_1d = data->d_collo_grad_1d;
158       }
159       void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v};
160       if (dim == 1) {
161         CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64;
162         elems_per_block         = elems_per_block > 0 ? elems_per_block : 1;
163         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
164         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
165         if (t_mode == CEED_TRANSPOSE) {
166           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
167         } else {
168           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
169         }
170       } else if (dim == 2) {
171         // Check if required threads is small enough to do multiple elems
172         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
173         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
174         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
175         if (t_mode == CEED_TRANSPOSE) {
176           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
177         } else {
178           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
179         }
180       } else if (dim == 3) {
181         const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1);
182         CeedInt       grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
183         CeedInt       shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
184         if (t_mode == CEED_TRANSPOSE) {
185           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
186         } else {
187           CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
188         }
189       }
190     } break;
191     case CEED_EVAL_WEIGHT: {
192       CeedInt Q_1d;
193       CeedInt block_size = data->block_sizes[2];
194       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
195       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
196       if (dim == 1) {
197         const CeedInt opt_elems       = block_size / Q_1d;
198         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
199         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
200         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
201       } else if (dim == 2) {
202         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
203         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
204         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
205         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
206       } else if (dim == 3) {
207         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
208         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
209         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
210         CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
211       }
212     } break;
213     // LCOV_EXCL_START
214     // Evaluate the divergence to/from the quadrature points
215     case CEED_EVAL_DIV:
216       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported");
217     // Evaluate the curl to/from the quadrature points
218     case CEED_EVAL_CURL:
219       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported");
220     // Take no action, BasisApply should not have been called
221     case CEED_EVAL_NONE:
222       return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context");
223       // LCOV_EXCL_STOP
224   }
225 
226   // Restore vectors
227   if (eval_mode != CEED_EVAL_WEIGHT) {
228     CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
229   }
230   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
231   return CEED_ERROR_SUCCESS;
232 }
233 
234 //------------------------------------------------------------------------------
235 // Destroy basis
236 //------------------------------------------------------------------------------
237 static int CeedBasisDestroy_Hip_shared(CeedBasis basis) {
238   Ceed ceed;
239   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
240 
241   CeedBasis_Hip_shared *data;
242   CeedCallBackend(CeedBasisGetData(basis, &data));
243 
244   CeedCallHip(ceed, hipModuleUnload(data->module));
245 
246   CeedCallHip(ceed, hipFree(data->d_q_weight_1d));
247   CeedCallHip(ceed, hipFree(data->d_interp_1d));
248   CeedCallHip(ceed, hipFree(data->d_grad_1d));
249   CeedCallHip(ceed, hipFree(data->d_collo_grad_1d));
250   CeedCallBackend(CeedFree(&data));
251 
252   return CEED_ERROR_SUCCESS;
253 }
254 
255 //------------------------------------------------------------------------------
256 // Create tensor basis
257 //------------------------------------------------------------------------------
258 int CeedBasisCreateTensorH1_Hip_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
259                                        const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
260   Ceed ceed;
261   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
262   CeedBasis_Hip_shared *data;
263   CeedCallBackend(CeedCalloc(1, &data));
264 
265   // Copy basis data to GPU
266   const CeedInt qBytes = Q_1d * sizeof(CeedScalar);
267   CeedCallHip(ceed, hipMalloc((void **)&data->d_q_weight_1d, qBytes));
268   CeedCallHip(ceed, hipMemcpy(data->d_q_weight_1d, q_weight_1d, qBytes, hipMemcpyHostToDevice));
269 
270   const CeedInt iBytes = qBytes * P_1d;
271   CeedCallHip(ceed, hipMalloc((void **)&data->d_interp_1d, iBytes));
272   CeedCallHip(ceed, hipMemcpy(data->d_interp_1d, interp_1d, iBytes, hipMemcpyHostToDevice));
273 
274   CeedCallHip(ceed, hipMalloc((void **)&data->d_grad_1d, iBytes));
275   CeedCallHip(ceed, hipMemcpy(data->d_grad_1d, grad_1d, iBytes, hipMemcpyHostToDevice));
276 
277   // Compute collocated gradient and copy to GPU
278   data->d_collo_grad_1d    = NULL;
279   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
280   if (has_collocated_grad) {
281     CeedScalar *collo_grad_1d;
282     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
283     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
284     CeedCallHip(ceed, hipMalloc((void **)&data->d_collo_grad_1d, qBytes * Q_1d));
285     CeedCallHip(ceed, hipMemcpy(data->d_collo_grad_1d, collo_grad_1d, qBytes * Q_1d, hipMemcpyHostToDevice));
286     CeedCallBackend(CeedFree(&collo_grad_1d));
287   }
288 
289   // Set number of threads per block for basis kernels
290   CeedInt num_comp;
291   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
292   CeedCallBackend(ComputeBasisThreadBlockSizes(dim, P_1d, Q_1d, num_comp, data->block_sizes));
293 
294   // Compile basis kernels
295   char *basis_kernel_path, *basis_kernel_source;
296   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor.h", &basis_kernel_path));
297   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
298   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
299   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n");
300   CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->module, 11, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
301                                   CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
302                                   "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_INTERP_BLOCK_SIZE", data->block_sizes[0], "BASIS_GRAD_BLOCK_SIZE",
303                                   data->block_sizes[1], "BASIS_WEIGHT_BLOCK_SIZE", data->block_sizes[2], "BASIS_HAS_COLLOCATED_GRAD",
304                                   has_collocated_grad));
305   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Interp", &data->Interp));
306   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
307   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Grad", &data->Grad));
308   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTranspose", &data->GradTranspose));
309   CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Weight", &data->Weight));
310   CeedCallBackend(CeedFree(&basis_kernel_path));
311   CeedCallBackend(CeedFree(&basis_kernel_source));
312 
313   CeedCallBackend(CeedBasisSetData(basis, data));
314 
315   // Register backend functions
316   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Hip_shared));
317   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Hip_shared));
318   return CEED_ERROR_SUCCESS;
319 }
320 
321 //------------------------------------------------------------------------------
322