1ff1e7120SSebastian Grimberg // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2ff1e7120SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ff1e7120SSebastian Grimberg // 4ff1e7120SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause 5ff1e7120SSebastian Grimberg // 6ff1e7120SSebastian Grimberg // This file is part of CEED: http://github.com/ceed 7ff1e7120SSebastian Grimberg 8ff1e7120SSebastian Grimberg #include <ceed.h> 9ff1e7120SSebastian Grimberg #include <ceed/backend.h> 10ff1e7120SSebastian Grimberg #include <ceed/jit-tools.h> 11ff1e7120SSebastian Grimberg #include <cuda.h> 12ff1e7120SSebastian Grimberg #include <cuda_runtime.h> 13ff1e7120SSebastian Grimberg #include <stdbool.h> 14ff1e7120SSebastian Grimberg #include <stddef.h> 15ff1e7120SSebastian Grimberg #include <string.h> 16ff1e7120SSebastian Grimberg 17ff1e7120SSebastian Grimberg #include "../cuda/ceed-cuda-common.h" 18ff1e7120SSebastian Grimberg #include "../cuda/ceed-cuda-compile.h" 19ff1e7120SSebastian Grimberg #include "ceed-cuda-ref.h" 20ff1e7120SSebastian Grimberg 21ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 22ff1e7120SSebastian Grimberg // Apply restriction 23ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 24ff1e7120SSebastian Grimberg static int CeedElemRestrictionApply_Cuda(CeedElemRestriction r, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) { 25ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 26ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(r, &impl)); 27ff1e7120SSebastian Grimberg Ceed ceed; 28ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed)); 29ff1e7120SSebastian Grimberg Ceed_Cuda *data; 30ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetData(ceed, &data)); 31ff1e7120SSebastian Grimberg const CeedInt warp_size = 32; 32ff1e7120SSebastian Grimberg const CeedInt block_size = warp_size; 33ff1e7120SSebastian Grimberg const CeedInt num_nodes = impl->num_nodes; 34ff1e7120SSebastian Grimberg CeedInt num_elem, elem_size; 35ff1e7120SSebastian Grimberg CeedElemRestrictionGetNumElements(r, &num_elem); 36ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size)); 37ff1e7120SSebastian Grimberg CUfunction kernel; 38ff1e7120SSebastian Grimberg 39ff1e7120SSebastian Grimberg // Get vectors 40ff1e7120SSebastian Grimberg const CeedScalar *d_u; 41ff1e7120SSebastian Grimberg CeedScalar *d_v; 42ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 43ff1e7120SSebastian Grimberg if (t_mode == CEED_TRANSPOSE) { 44ff1e7120SSebastian Grimberg // Sum into for transpose mode, e-vec to l-vec 45ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 46ff1e7120SSebastian Grimberg } else { 47ff1e7120SSebastian Grimberg // Overwrite for notranspose mode, l-vec to e-vec 48ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 49ff1e7120SSebastian Grimberg } 50ff1e7120SSebastian Grimberg 51ff1e7120SSebastian Grimberg // Restrict 52ff1e7120SSebastian Grimberg if (t_mode == CEED_NOTRANSPOSE) { 53ff1e7120SSebastian Grimberg // L-vector -> E-vector 54ff1e7120SSebastian Grimberg if (impl->d_ind) { 55ff1e7120SSebastian Grimberg // -- Offsets provided 56ff1e7120SSebastian Grimberg kernel = impl->OffsetNoTranspose; 57ff1e7120SSebastian Grimberg void *args[] = {&num_elem, &impl->d_ind, &d_u, &d_v}; 58ff1e7120SSebastian Grimberg CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024; 59ff1e7120SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args)); 60ff1e7120SSebastian Grimberg } else { 61ff1e7120SSebastian Grimberg // -- Strided restriction 62ff1e7120SSebastian Grimberg kernel = impl->StridedNoTranspose; 63ff1e7120SSebastian Grimberg void *args[] = {&num_elem, &d_u, &d_v}; 64ff1e7120SSebastian Grimberg CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024; 65ff1e7120SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args)); 66ff1e7120SSebastian Grimberg } 67ff1e7120SSebastian Grimberg } else { 68ff1e7120SSebastian Grimberg // E-vector -> L-vector 69ff1e7120SSebastian Grimberg if (impl->d_ind) { 70ff1e7120SSebastian Grimberg // -- Offsets provided 71ff1e7120SSebastian Grimberg kernel = impl->OffsetTranspose; 72ff1e7120SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 73ff1e7120SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args)); 74ff1e7120SSebastian Grimberg } else { 75ff1e7120SSebastian Grimberg // -- Strided restriction 76ff1e7120SSebastian Grimberg kernel = impl->StridedTranspose; 77ff1e7120SSebastian Grimberg void *args[] = {&num_elem, &d_u, &d_v}; 78ff1e7120SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, kernel, CeedDivUpInt(num_nodes, block_size), block_size, args)); 79ff1e7120SSebastian Grimberg } 80ff1e7120SSebastian Grimberg } 81ff1e7120SSebastian Grimberg 82ff1e7120SSebastian Grimberg if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL; 83ff1e7120SSebastian Grimberg 84ff1e7120SSebastian Grimberg // Restore arrays 85ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 86ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 87ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 88ff1e7120SSebastian Grimberg } 89ff1e7120SSebastian Grimberg 90ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 91ff1e7120SSebastian Grimberg // Get offsets 92ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 93ff1e7120SSebastian Grimberg static int CeedElemRestrictionGetOffsets_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) { 94ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 95ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 96ff1e7120SSebastian Grimberg 97ff1e7120SSebastian Grimberg switch (mem_type) { 98ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 99ff1e7120SSebastian Grimberg *offsets = impl->h_ind; 100ff1e7120SSebastian Grimberg break; 101ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 102ff1e7120SSebastian Grimberg *offsets = impl->d_ind; 103ff1e7120SSebastian Grimberg break; 104ff1e7120SSebastian Grimberg } 105ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 106ff1e7120SSebastian Grimberg } 107ff1e7120SSebastian Grimberg 108ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 109ff1e7120SSebastian Grimberg // Destroy restriction 110ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 111ff1e7120SSebastian Grimberg static int CeedElemRestrictionDestroy_Cuda(CeedElemRestriction r) { 112ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 113ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(r, &impl)); 114ff1e7120SSebastian Grimberg 115ff1e7120SSebastian Grimberg Ceed ceed; 116ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed)); 117ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cuModuleUnload(impl->module)); 118ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_ind_allocated)); 119ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_ind_allocated)); 120ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_t_offsets)); 121ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_t_indices)); 122ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_l_vec_indices)); 123ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 124ff1e7120SSebastian Grimberg 125ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 126ff1e7120SSebastian Grimberg } 127ff1e7120SSebastian Grimberg 128ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 129ff1e7120SSebastian Grimberg // Create transpose offsets and indices 130ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 131ff1e7120SSebastian Grimberg static int CeedElemRestrictionOffset_Cuda(const CeedElemRestriction r, const CeedInt *indices) { 132ff1e7120SSebastian Grimberg Ceed ceed; 133ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed)); 134ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 135ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(r, &impl)); 136ff1e7120SSebastian Grimberg CeedSize l_size; 137ff1e7120SSebastian Grimberg CeedInt num_elem, elem_size, num_comp; 138ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem)); 139ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size)); 140ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetLVectorSize(r, &l_size)); 141ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp)); 142ff1e7120SSebastian Grimberg 143ff1e7120SSebastian Grimberg // Count num_nodes 144ff1e7120SSebastian Grimberg bool *is_node; 145ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &is_node)); 146ff1e7120SSebastian Grimberg const CeedInt size_indices = num_elem * elem_size; 147ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1; 148ff1e7120SSebastian Grimberg CeedInt num_nodes = 0; 149ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i]; 150ff1e7120SSebastian Grimberg impl->num_nodes = num_nodes; 151ff1e7120SSebastian Grimberg 152ff1e7120SSebastian Grimberg // L-vector offsets array 153ff1e7120SSebastian Grimberg CeedInt *ind_to_offset, *l_vec_indices; 154ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &ind_to_offset)); 155ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices)); 156ff1e7120SSebastian Grimberg CeedInt j = 0; 157ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < l_size; i++) { 158ff1e7120SSebastian Grimberg if (is_node[i]) { 159ff1e7120SSebastian Grimberg l_vec_indices[j] = i; 160ff1e7120SSebastian Grimberg ind_to_offset[i] = j++; 161ff1e7120SSebastian Grimberg } 162ff1e7120SSebastian Grimberg } 163ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&is_node)); 164ff1e7120SSebastian Grimberg 165ff1e7120SSebastian Grimberg // Compute transpose offsets and indices 166ff1e7120SSebastian Grimberg const CeedInt size_offsets = num_nodes + 1; 167ff1e7120SSebastian Grimberg CeedInt *t_offsets; 168ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(size_offsets, &t_offsets)); 169ff1e7120SSebastian Grimberg CeedInt *t_indices; 170ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(size_indices, &t_indices)); 171ff1e7120SSebastian Grimberg // Count node multiplicity 172ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 173ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1]; 174ff1e7120SSebastian Grimberg } 175ff1e7120SSebastian Grimberg // Convert to running sum 176ff1e7120SSebastian Grimberg for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1]; 177ff1e7120SSebastian Grimberg // List all E-vec indices associated with L-vec node 178ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 179ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) { 180ff1e7120SSebastian Grimberg const CeedInt lid = elem_size * e + i; 181ff1e7120SSebastian Grimberg const CeedInt gid = indices[lid]; 182ff1e7120SSebastian Grimberg t_indices[t_offsets[ind_to_offset[gid]]++] = lid; 183ff1e7120SSebastian Grimberg } 184ff1e7120SSebastian Grimberg } 185ff1e7120SSebastian Grimberg // Reset running sum 186ff1e7120SSebastian Grimberg for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1]; 187ff1e7120SSebastian Grimberg t_offsets[0] = 0; 188ff1e7120SSebastian Grimberg 189ff1e7120SSebastian Grimberg // Copy data to device 190ff1e7120SSebastian Grimberg // -- L-vector indices 191ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt))); 192ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), cudaMemcpyHostToDevice)); 193ff1e7120SSebastian Grimberg // -- Transpose offsets 194ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt))); 195ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), cudaMemcpyHostToDevice)); 196ff1e7120SSebastian Grimberg // -- Transpose indices 197ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt))); 198ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), cudaMemcpyHostToDevice)); 199ff1e7120SSebastian Grimberg 200ff1e7120SSebastian Grimberg // Cleanup 201ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&ind_to_offset)); 202ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&l_vec_indices)); 203ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_offsets)); 204ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_indices)); 205ff1e7120SSebastian Grimberg 206ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 207ff1e7120SSebastian Grimberg } 208ff1e7120SSebastian Grimberg 209ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 210ff1e7120SSebastian Grimberg // Create restriction 211ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 212fcbe8c06SSebastian Grimberg int CeedElemRestrictionCreate_Cuda(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *indices, const bool *orients, 213*0c73c039SSebastian Grimberg const CeedInt8 *curl_orients, CeedElemRestriction r) { 214ff1e7120SSebastian Grimberg Ceed ceed; 215ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(r, &ceed)); 216ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 217ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 218ff1e7120SSebastian Grimberg CeedInt num_elem, num_comp, elem_size; 219ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(r, &num_elem)); 220ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(r, &num_comp)); 221ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(r, &elem_size)); 222ff1e7120SSebastian Grimberg CeedInt size = num_elem * elem_size; 223ff1e7120SSebastian Grimberg CeedInt strides[3] = {1, size, elem_size}; 224ff1e7120SSebastian Grimberg CeedInt comp_stride = 1; 225ff1e7120SSebastian Grimberg 226ff1e7120SSebastian Grimberg // Stride data 227ff1e7120SSebastian Grimberg bool is_strided; 228ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionIsStrided(r, &is_strided)); 229ff1e7120SSebastian Grimberg if (is_strided) { 230ff1e7120SSebastian Grimberg bool has_backend_strides; 231ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionHasBackendStrides(r, &has_backend_strides)); 232ff1e7120SSebastian Grimberg if (!has_backend_strides) { 233ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetStrides(r, &strides)); 234ff1e7120SSebastian Grimberg } 235ff1e7120SSebastian Grimberg } else { 236ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCompStride(r, &comp_stride)); 237ff1e7120SSebastian Grimberg } 238ff1e7120SSebastian Grimberg 239ff1e7120SSebastian Grimberg impl->h_ind = NULL; 240ff1e7120SSebastian Grimberg impl->h_ind_allocated = NULL; 241ff1e7120SSebastian Grimberg impl->d_ind = NULL; 242ff1e7120SSebastian Grimberg impl->d_ind_allocated = NULL; 243ff1e7120SSebastian Grimberg impl->d_t_indices = NULL; 244ff1e7120SSebastian Grimberg impl->d_t_offsets = NULL; 245ff1e7120SSebastian Grimberg impl->num_nodes = size; 246ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetData(r, impl)); 247ff1e7120SSebastian Grimberg CeedInt layout[3] = {1, elem_size * num_elem, elem_size}; 248ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetELayout(r, layout)); 249ff1e7120SSebastian Grimberg 250ff1e7120SSebastian Grimberg // Set up device indices/offset arrays 251ff1e7120SSebastian Grimberg switch (mem_type) { 252ff1e7120SSebastian Grimberg case CEED_MEM_HOST: { 253ff1e7120SSebastian Grimberg switch (copy_mode) { 254ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 255ff1e7120SSebastian Grimberg impl->h_ind_allocated = (CeedInt *)indices; 256ff1e7120SSebastian Grimberg impl->h_ind = (CeedInt *)indices; 257ff1e7120SSebastian Grimberg break; 258ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 259ff1e7120SSebastian Grimberg impl->h_ind = (CeedInt *)indices; 260ff1e7120SSebastian Grimberg break; 261ff1e7120SSebastian Grimberg case CEED_COPY_VALUES: 262ff1e7120SSebastian Grimberg if (indices != NULL) { 263ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated)); 264ff1e7120SSebastian Grimberg memcpy(impl->h_ind_allocated, indices, elem_size * num_elem * sizeof(CeedInt)); 265ff1e7120SSebastian Grimberg impl->h_ind = impl->h_ind_allocated; 266ff1e7120SSebastian Grimberg } 267ff1e7120SSebastian Grimberg break; 268ff1e7120SSebastian Grimberg } 269ff1e7120SSebastian Grimberg if (indices != NULL) { 270ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_ind, size * sizeof(CeedInt))); 271ff1e7120SSebastian Grimberg impl->d_ind_allocated = impl->d_ind; // We own the device memory 272ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), cudaMemcpyHostToDevice)); 273ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionOffset_Cuda(r, indices)); 274ff1e7120SSebastian Grimberg } 275ff1e7120SSebastian Grimberg break; 276ff1e7120SSebastian Grimberg } 277ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: { 278ff1e7120SSebastian Grimberg switch (copy_mode) { 279ff1e7120SSebastian Grimberg case CEED_COPY_VALUES: 280ff1e7120SSebastian Grimberg if (indices != NULL) { 281ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_ind, size * sizeof(CeedInt))); 282ff1e7120SSebastian Grimberg impl->d_ind_allocated = impl->d_ind; // We own the device memory 283ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_ind, indices, size * sizeof(CeedInt), cudaMemcpyDeviceToDevice)); 284ff1e7120SSebastian Grimberg } 285ff1e7120SSebastian Grimberg break; 286ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 287ff1e7120SSebastian Grimberg impl->d_ind = (CeedInt *)indices; 288ff1e7120SSebastian Grimberg impl->d_ind_allocated = impl->d_ind; 289ff1e7120SSebastian Grimberg break; 290ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 291ff1e7120SSebastian Grimberg impl->d_ind = (CeedInt *)indices; 292ff1e7120SSebastian Grimberg } 293ff1e7120SSebastian Grimberg if (indices != NULL) { 294ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(elem_size * num_elem, &impl->h_ind_allocated)); 295ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->h_ind_allocated, impl->d_ind, elem_size * num_elem * sizeof(CeedInt), cudaMemcpyDeviceToHost)); 296ff1e7120SSebastian Grimberg impl->h_ind = impl->h_ind_allocated; 297ff1e7120SSebastian Grimberg CeedCallBackend(CeedElemRestrictionOffset_Cuda(r, indices)); 298ff1e7120SSebastian Grimberg } 299ff1e7120SSebastian Grimberg break; 300ff1e7120SSebastian Grimberg } 301ff1e7120SSebastian Grimberg // LCOV_EXCL_START 302ff1e7120SSebastian Grimberg default: 303ff1e7120SSebastian Grimberg return CeedError(ceed, CEED_ERROR_BACKEND, "Only MemType = HOST or DEVICE supported"); 304ff1e7120SSebastian Grimberg // LCOV_EXCL_STOP 305ff1e7120SSebastian Grimberg } 306ff1e7120SSebastian Grimberg 307ff1e7120SSebastian Grimberg // Compile CUDA kernels 308ff1e7120SSebastian Grimberg CeedInt num_nodes = impl->num_nodes; 309ff1e7120SSebastian Grimberg char *restriction_kernel_path, *restriction_kernel_source; 310ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction.h", &restriction_kernel_path)); 31123d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 312ff1e7120SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 31323d4529eSJeremy L Thompson CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 314ff1e7120SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 8, "RESTR_ELEM_SIZE", elem_size, "RESTR_NUM_ELEM", num_elem, 315ff1e7120SSebastian Grimberg "RESTR_NUM_COMP", num_comp, "RESTR_NUM_NODES", num_nodes, "RESTR_COMP_STRIDE", comp_stride, "RESTR_STRIDE_NODES", 316ff1e7120SSebastian Grimberg strides[0], "RESTR_STRIDE_COMP", strides[1], "RESTR_STRIDE_ELEM", strides[2])); 317ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedTranspose", &impl->StridedTranspose)); 318ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedNoTranspose", &impl->StridedNoTranspose)); 319ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->OffsetTranspose)); 320ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->OffsetNoTranspose)); 321ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_path)); 322ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_source)); 323ff1e7120SSebastian Grimberg 324ff1e7120SSebastian Grimberg // Register backend functions 325ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Apply", CeedElemRestrictionApply_Cuda)); 326ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "ApplyUnsigned", CeedElemRestrictionApply_Cuda)); 327ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "GetOffsets", CeedElemRestrictionGetOffsets_Cuda)); 328ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", r, "Destroy", CeedElemRestrictionDestroy_Cuda)); 329ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 330ff1e7120SSebastian Grimberg } 331ff1e7120SSebastian Grimberg 332ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 333