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 //------------------------------------------------------------------------------ 22cf8cbdd6SSebastian Grimberg // Compile restriction kernels 23cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 24cf8cbdd6SSebastian Grimberg static inline int CeedElemRestrictionSetupCompile_Cuda(CeedElemRestriction rstr) { 25cf8cbdd6SSebastian Grimberg Ceed ceed; 26cf8cbdd6SSebastian Grimberg bool is_deterministic; 2722070f95SJeremy L Thompson char *restriction_kernel_source; 2822070f95SJeremy L Thompson const char *restriction_kernel_path; 29cf8cbdd6SSebastian Grimberg CeedInt num_elem, num_comp, elem_size, comp_stride; 30cf8cbdd6SSebastian Grimberg CeedRestrictionType rstr_type; 31cf8cbdd6SSebastian Grimberg CeedElemRestriction_Cuda *impl; 32cf8cbdd6SSebastian Grimberg 33cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 34cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 35cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 36cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 37cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 38cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride)); 39cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 40cf8cbdd6SSebastian Grimberg is_deterministic = impl->d_l_vec_indices != NULL; 41cf8cbdd6SSebastian Grimberg 42cf8cbdd6SSebastian Grimberg // Compile CUDA kernels 43cf8cbdd6SSebastian Grimberg switch (rstr_type) { 44cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 45cf8cbdd6SSebastian Grimberg CeedInt strides[3] = {1, num_elem * elem_size, elem_size}; 46cf8cbdd6SSebastian Grimberg bool has_backend_strides; 47cf8cbdd6SSebastian Grimberg 48cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides)); 49cf8cbdd6SSebastian Grimberg if (!has_backend_strides) { 5056c48462SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(rstr, strides)); 51cf8cbdd6SSebastian Grimberg } 52cf8cbdd6SSebastian Grimberg 53cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-strided.h", &restriction_kernel_path)); 54cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 55cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 56cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 57cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 58cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_STRIDE_NODES", strides[0], "RSTR_STRIDE_COMP", strides[1], "RSTR_STRIDE_ELEM", 59cf8cbdd6SSebastian Grimberg strides[2])); 60cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedNoTranspose", &impl->ApplyNoTranspose)); 61cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedTranspose", &impl->ApplyTranspose)); 62cf8cbdd6SSebastian Grimberg } break; 63cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 64cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &restriction_kernel_path)); 65cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 66cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 67cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 68cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 69cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 70cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 71cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyNoTranspose)); 72cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyTranspose)); 73cf8cbdd6SSebastian Grimberg } break; 74cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 7522070f95SJeremy L Thompson const char *offset_kernel_path; 76cf8cbdd6SSebastian Grimberg 77cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-oriented.h", &restriction_kernel_path)); 78cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 79cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 80cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &offset_kernel_path)); 81cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToInitializedBuffer(ceed, offset_kernel_path, &restriction_kernel_source)); 82cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 83cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 84cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 85cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 86cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedNoTranspose", &impl->ApplyNoTranspose)); 87cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnsignedNoTranspose)); 88cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedTranspose", &impl->ApplyTranspose)); 89cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnsignedTranspose)); 90cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&offset_kernel_path)); 91cf8cbdd6SSebastian Grimberg } break; 92cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 9322070f95SJeremy L Thompson const char *offset_kernel_path; 94cf8cbdd6SSebastian Grimberg 95cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-curl-oriented.h", &restriction_kernel_path)); 96cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 97cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 98cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &offset_kernel_path)); 99cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToInitializedBuffer(ceed, offset_kernel_path, &restriction_kernel_source)); 100cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 101cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 102cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 103cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 104cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedNoTranspose", &impl->ApplyNoTranspose)); 105cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedNoTranspose", &impl->ApplyUnsignedNoTranspose)); 106cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnorientedNoTranspose)); 107cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedTranspose", &impl->ApplyTranspose)); 108cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedTranspose", &impl->ApplyUnsignedTranspose)); 109cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnorientedTranspose)); 110cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&offset_kernel_path)); 111cf8cbdd6SSebastian Grimberg } break; 112cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_POINTS: { 113cf8cbdd6SSebastian Grimberg // LCOV_EXCL_START 114cf8cbdd6SSebastian Grimberg return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement restriction CeedElemRestrictionAtPoints"); 115cf8cbdd6SSebastian Grimberg // LCOV_EXCL_STOP 116cf8cbdd6SSebastian Grimberg } break; 117cf8cbdd6SSebastian Grimberg } 118cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_path)); 119cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_source)); 120cf8cbdd6SSebastian Grimberg return CEED_ERROR_SUCCESS; 121cf8cbdd6SSebastian Grimberg } 122cf8cbdd6SSebastian Grimberg 123cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 124dce49693SSebastian Grimberg // Core apply restriction code 125ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 126dce49693SSebastian Grimberg static inline int CeedElemRestrictionApply_Cuda_Core(CeedElemRestriction rstr, CeedTransposeMode t_mode, bool use_signs, bool use_orients, 127dce49693SSebastian Grimberg CeedVector u, CeedVector v, CeedRequest *request) { 128ff1e7120SSebastian Grimberg Ceed ceed; 129dce49693SSebastian Grimberg CeedRestrictionType rstr_type; 130ff1e7120SSebastian Grimberg const CeedScalar *d_u; 131ff1e7120SSebastian Grimberg CeedScalar *d_v; 132ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 133ca735530SJeremy L Thompson 134dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 135dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 136dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 137cf8cbdd6SSebastian Grimberg 138cf8cbdd6SSebastian Grimberg // Assemble kernel if needed 139cf8cbdd6SSebastian Grimberg if (!impl->module) { 140cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetupCompile_Cuda(rstr)); 141cf8cbdd6SSebastian Grimberg } 142ca735530SJeremy L Thompson 143ca735530SJeremy L Thompson // Get vectors 144ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 145ff1e7120SSebastian Grimberg if (t_mode == CEED_TRANSPOSE) { 146ff1e7120SSebastian Grimberg // Sum into for transpose mode, e-vec to l-vec 147ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 148ff1e7120SSebastian Grimberg } else { 149ff1e7120SSebastian Grimberg // Overwrite for notranspose mode, l-vec to e-vec 150ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 151ff1e7120SSebastian Grimberg } 152ff1e7120SSebastian Grimberg 153ff1e7120SSebastian Grimberg // Restrict 154ff1e7120SSebastian Grimberg if (t_mode == CEED_NOTRANSPOSE) { 155ff1e7120SSebastian Grimberg // L-vector -> E-vector 156cf8cbdd6SSebastian Grimberg CeedInt elem_size; 157cf8cbdd6SSebastian Grimberg 158cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 159dce49693SSebastian Grimberg const CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024; 160cf8cbdd6SSebastian Grimberg const CeedInt grid = CeedDivUpInt(impl->num_nodes, block_size); 16158549094SSebastian Grimberg 162dce49693SSebastian Grimberg switch (rstr_type) { 163dce49693SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 164cf8cbdd6SSebastian Grimberg void *args[] = {&d_u, &d_v}; 16558549094SSebastian Grimberg 166cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 167dce49693SSebastian Grimberg } break; 168dce49693SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 169*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 170dce49693SSebastian Grimberg 171cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 172dce49693SSebastian Grimberg } break; 173dce49693SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 174dce49693SSebastian Grimberg if (use_signs) { 175*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v}; 176dce49693SSebastian Grimberg 177cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 178dce49693SSebastian Grimberg } else { 179*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 180dce49693SSebastian Grimberg 181cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args)); 182dce49693SSebastian Grimberg } 183dce49693SSebastian Grimberg } break; 184dce49693SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 185dce49693SSebastian Grimberg if (use_signs && use_orients) { 186*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 187dce49693SSebastian Grimberg 188cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 189dce49693SSebastian Grimberg } else if (use_orients) { 190*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 191dce49693SSebastian Grimberg 192cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args)); 193dce49693SSebastian Grimberg } else { 194*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 195dce49693SSebastian Grimberg 196cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedNoTranspose, grid, block_size, args)); 197dce49693SSebastian Grimberg } 198dce49693SSebastian Grimberg } break; 199b3d03e38SSebastian Grimberg case CEED_RESTRICTION_POINTS: { 200b3d03e38SSebastian Grimberg // LCOV_EXCL_START 201b3d03e38SSebastian Grimberg return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement restriction CeedElemRestrictionAtPoints"); 202b3d03e38SSebastian Grimberg // LCOV_EXCL_STOP 203b3d03e38SSebastian Grimberg } break; 204ff1e7120SSebastian Grimberg } 205ff1e7120SSebastian Grimberg } else { 206ff1e7120SSebastian Grimberg // E-vector -> L-vector 207cf8cbdd6SSebastian Grimberg const bool is_deterministic = impl->d_l_vec_indices != NULL; 208dce49693SSebastian Grimberg const CeedInt block_size = 32; 209cf8cbdd6SSebastian Grimberg const CeedInt grid = CeedDivUpInt(impl->num_nodes, block_size); 210ca735530SJeremy L Thompson 211dce49693SSebastian Grimberg switch (rstr_type) { 212dce49693SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 213cf8cbdd6SSebastian Grimberg void *args[] = {&d_u, &d_v}; 214dce49693SSebastian Grimberg 215cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 216dce49693SSebastian Grimberg } break; 217dce49693SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 218cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 219*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 22058549094SSebastian Grimberg 221cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 222ff1e7120SSebastian Grimberg } else { 22358549094SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 22458549094SSebastian Grimberg 225cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 22658549094SSebastian Grimberg } 227dce49693SSebastian Grimberg } break; 228dce49693SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 229dce49693SSebastian Grimberg if (use_signs) { 230cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 231*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v}; 23258549094SSebastian Grimberg 233cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 234dce49693SSebastian Grimberg } else { 2357aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_orients, &d_u, &d_v}; 2367aa91133SSebastian Grimberg 237cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2387aa91133SSebastian Grimberg } 2397aa91133SSebastian Grimberg } else { 240cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 241*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 242dce49693SSebastian Grimberg 243cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 244dce49693SSebastian Grimberg } else { 245dce49693SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 246dce49693SSebastian Grimberg 247cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 248dce49693SSebastian Grimberg } 249dce49693SSebastian Grimberg } 250dce49693SSebastian Grimberg } break; 251dce49693SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 252dce49693SSebastian Grimberg if (use_signs && use_orients) { 253cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 254*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 255dce49693SSebastian Grimberg 256cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2577aa91133SSebastian Grimberg } else { 2587aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v}; 2597aa91133SSebastian Grimberg 260cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2617aa91133SSebastian Grimberg } 262dce49693SSebastian Grimberg } else if (use_orients) { 263cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 264*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 265dce49693SSebastian Grimberg 266cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 267dce49693SSebastian Grimberg } else { 2687aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v}; 2697aa91133SSebastian Grimberg 270cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 2717aa91133SSebastian Grimberg } 2727aa91133SSebastian Grimberg } else { 273cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 274*a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 275dce49693SSebastian Grimberg 276cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args)); 277dce49693SSebastian Grimberg } else { 278dce49693SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 279dce49693SSebastian Grimberg 280cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args)); 281dce49693SSebastian Grimberg } 282dce49693SSebastian Grimberg } 283dce49693SSebastian Grimberg } break; 284b3d03e38SSebastian Grimberg case CEED_RESTRICTION_POINTS: { 285b3d03e38SSebastian Grimberg // LCOV_EXCL_START 286b3d03e38SSebastian Grimberg return CeedError(ceed, CEED_ERROR_UNSUPPORTED, "Backend does not implement restriction CeedElemRestrictionAtPoints"); 287b3d03e38SSebastian Grimberg // LCOV_EXCL_STOP 288b3d03e38SSebastian Grimberg } break; 289ff1e7120SSebastian Grimberg } 290ff1e7120SSebastian Grimberg } 291ff1e7120SSebastian Grimberg 292ff1e7120SSebastian Grimberg if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL; 293ff1e7120SSebastian Grimberg 294ff1e7120SSebastian Grimberg // Restore arrays 295ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 296ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 297ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 298ff1e7120SSebastian Grimberg } 299ff1e7120SSebastian Grimberg 300ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 301dce49693SSebastian Grimberg // Apply restriction 302dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 303dce49693SSebastian Grimberg static int CeedElemRestrictionApply_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) { 304dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, true, true, u, v, request); 305dce49693SSebastian Grimberg } 306dce49693SSebastian Grimberg 307dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 308dce49693SSebastian Grimberg // Apply unsigned restriction 309dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 310dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnsigned_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, 311dce49693SSebastian Grimberg CeedRequest *request) { 312dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, true, u, v, request); 313dce49693SSebastian Grimberg } 314dce49693SSebastian Grimberg 315dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 316dce49693SSebastian Grimberg // Apply unoriented restriction 317dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 318dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnoriented_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, 319dce49693SSebastian Grimberg CeedRequest *request) { 320dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, false, u, v, request); 321dce49693SSebastian Grimberg } 322dce49693SSebastian Grimberg 323dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 324ff1e7120SSebastian Grimberg // Get offsets 325ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 326ff1e7120SSebastian Grimberg static int CeedElemRestrictionGetOffsets_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) { 327ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 328ff1e7120SSebastian Grimberg 329ca735530SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 330ff1e7120SSebastian Grimberg switch (mem_type) { 331ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 332*a267acd1SJeremy L Thompson *offsets = impl->h_offsets; 333ff1e7120SSebastian Grimberg break; 334ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 335*a267acd1SJeremy L Thompson *offsets = impl->d_offsets; 336ff1e7120SSebastian Grimberg break; 337ff1e7120SSebastian Grimberg } 338ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 339ff1e7120SSebastian Grimberg } 340ff1e7120SSebastian Grimberg 341ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 342dce49693SSebastian Grimberg // Get orientations 343dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 344dce49693SSebastian Grimberg static int CeedElemRestrictionGetOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const bool **orients) { 345dce49693SSebastian Grimberg CeedElemRestriction_Cuda *impl; 346dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 347dce49693SSebastian Grimberg 348dce49693SSebastian Grimberg switch (mem_type) { 349dce49693SSebastian Grimberg case CEED_MEM_HOST: 350dce49693SSebastian Grimberg *orients = impl->h_orients; 351dce49693SSebastian Grimberg break; 352dce49693SSebastian Grimberg case CEED_MEM_DEVICE: 353dce49693SSebastian Grimberg *orients = impl->d_orients; 354dce49693SSebastian Grimberg break; 355dce49693SSebastian Grimberg } 356dce49693SSebastian Grimberg return CEED_ERROR_SUCCESS; 357dce49693SSebastian Grimberg } 358dce49693SSebastian Grimberg 359dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 360dce49693SSebastian Grimberg // Get curl-conforming orientations 361dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 362dce49693SSebastian Grimberg static int CeedElemRestrictionGetCurlOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt8 **curl_orients) { 363dce49693SSebastian Grimberg CeedElemRestriction_Cuda *impl; 364dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 365dce49693SSebastian Grimberg 366dce49693SSebastian Grimberg switch (mem_type) { 367dce49693SSebastian Grimberg case CEED_MEM_HOST: 368dce49693SSebastian Grimberg *curl_orients = impl->h_curl_orients; 369dce49693SSebastian Grimberg break; 370dce49693SSebastian Grimberg case CEED_MEM_DEVICE: 371dce49693SSebastian Grimberg *curl_orients = impl->d_curl_orients; 372dce49693SSebastian Grimberg break; 373dce49693SSebastian Grimberg } 374dce49693SSebastian Grimberg return CEED_ERROR_SUCCESS; 375dce49693SSebastian Grimberg } 376dce49693SSebastian Grimberg 377dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 378ff1e7120SSebastian Grimberg // Destroy restriction 379ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 380dce49693SSebastian Grimberg static int CeedElemRestrictionDestroy_Cuda(CeedElemRestriction rstr) { 381ff1e7120SSebastian Grimberg Ceed ceed; 382ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 383ca735530SJeremy L Thompson 384dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 385dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 386cf8cbdd6SSebastian Grimberg if (impl->module) { 387ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cuModuleUnload(impl->module)); 388cf8cbdd6SSebastian Grimberg } 389*a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_offsets_owned)); 390*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaFree(impl->d_offsets_owned)); 391ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_t_offsets)); 392ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_t_indices)); 393ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_l_vec_indices)); 394*a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_orients_owned)); 395*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaFree(impl->d_orients_owned)); 396*a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_curl_orients_owned)); 397*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaFree(impl->d_curl_orients_owned)); 398ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 399ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 400ff1e7120SSebastian Grimberg } 401ff1e7120SSebastian Grimberg 402ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 403ff1e7120SSebastian Grimberg // Create transpose offsets and indices 404ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 405dce49693SSebastian Grimberg static int CeedElemRestrictionOffset_Cuda(const CeedElemRestriction rstr, const CeedInt *indices) { 406ff1e7120SSebastian Grimberg Ceed ceed; 407ca735530SJeremy L Thompson bool *is_node; 408ff1e7120SSebastian Grimberg CeedSize l_size; 409ca735530SJeremy L Thompson CeedInt num_elem, elem_size, num_comp, num_nodes = 0; 410ca735530SJeremy L Thompson CeedInt *ind_to_offset, *l_vec_indices, *t_offsets, *t_indices; 411ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 412ca735530SJeremy L Thompson 413dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 414dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 415dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 416dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 417dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size)); 418dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 419ca735530SJeremy L Thompson const CeedInt size_indices = num_elem * elem_size; 420ff1e7120SSebastian Grimberg 421ff1e7120SSebastian Grimberg // Count num_nodes 422ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &is_node)); 423ca735530SJeremy L Thompson 424ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1; 425ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i]; 426ff1e7120SSebastian Grimberg impl->num_nodes = num_nodes; 427ff1e7120SSebastian Grimberg 428ff1e7120SSebastian Grimberg // L-vector offsets array 429ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &ind_to_offset)); 430ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices)); 431ca735530SJeremy L Thompson for (CeedInt i = 0, j = 0; i < l_size; i++) { 432ff1e7120SSebastian Grimberg if (is_node[i]) { 433ff1e7120SSebastian Grimberg l_vec_indices[j] = i; 434ff1e7120SSebastian Grimberg ind_to_offset[i] = j++; 435ff1e7120SSebastian Grimberg } 436ff1e7120SSebastian Grimberg } 437ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&is_node)); 438ff1e7120SSebastian Grimberg 439ff1e7120SSebastian Grimberg // Compute transpose offsets and indices 440ff1e7120SSebastian Grimberg const CeedInt size_offsets = num_nodes + 1; 441ca735530SJeremy L Thompson 442ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(size_offsets, &t_offsets)); 443ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(size_indices, &t_indices)); 444ff1e7120SSebastian Grimberg // Count node multiplicity 445ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 446ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1]; 447ff1e7120SSebastian Grimberg } 448ff1e7120SSebastian Grimberg // Convert to running sum 449ff1e7120SSebastian Grimberg for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1]; 450ff1e7120SSebastian Grimberg // List all E-vec indices associated with L-vec node 451ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 452ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) { 453ff1e7120SSebastian Grimberg const CeedInt lid = elem_size * e + i; 454ff1e7120SSebastian Grimberg const CeedInt gid = indices[lid]; 455ca735530SJeremy L Thompson 456ff1e7120SSebastian Grimberg t_indices[t_offsets[ind_to_offset[gid]]++] = lid; 457ff1e7120SSebastian Grimberg } 458ff1e7120SSebastian Grimberg } 459ff1e7120SSebastian Grimberg // Reset running sum 460ff1e7120SSebastian Grimberg for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1]; 461ff1e7120SSebastian Grimberg t_offsets[0] = 0; 462ff1e7120SSebastian Grimberg 463ff1e7120SSebastian Grimberg // Copy data to device 464ff1e7120SSebastian Grimberg // -- L-vector indices 465ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt))); 466ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), cudaMemcpyHostToDevice)); 467ff1e7120SSebastian Grimberg // -- Transpose offsets 468ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt))); 469ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), cudaMemcpyHostToDevice)); 470ff1e7120SSebastian Grimberg // -- Transpose indices 471ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt))); 472ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), cudaMemcpyHostToDevice)); 473ff1e7120SSebastian Grimberg 474ff1e7120SSebastian Grimberg // Cleanup 475ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&ind_to_offset)); 476ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&l_vec_indices)); 477ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_offsets)); 478ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_indices)); 479ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 480ff1e7120SSebastian Grimberg } 481ff1e7120SSebastian Grimberg 482ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 483ff1e7120SSebastian Grimberg // Create restriction 484ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 485*a267acd1SJeremy L Thompson int CeedElemRestrictionCreate_Cuda(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *offsets, const bool *orients, 486dce49693SSebastian Grimberg const CeedInt8 *curl_orients, CeedElemRestriction rstr) { 487ca735530SJeremy L Thompson Ceed ceed, ceed_parent; 488dce49693SSebastian Grimberg bool is_deterministic; 489cf8cbdd6SSebastian Grimberg CeedInt num_elem, elem_size; 490ca735530SJeremy L Thompson CeedRestrictionType rstr_type; 491ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 492ca735530SJeremy L Thompson 493dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 494ca735530SJeremy L Thompson CeedCallBackend(CeedGetParent(ceed, &ceed_parent)); 495ca735530SJeremy L Thompson CeedCallBackend(CeedIsDeterministic(ceed_parent, &is_deterministic)); 496dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 497dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 49822eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 499ca735530SJeremy L Thompson const CeedInt size = num_elem * elem_size; 500ff1e7120SSebastian Grimberg 501dce49693SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 502dce49693SSebastian Grimberg impl->num_nodes = size; 503dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetData(rstr, impl)); 50422eb1385SJeremy L Thompson 50522eb1385SJeremy L Thompson // Set layouts 50622eb1385SJeremy L Thompson { 50722eb1385SJeremy L Thompson bool has_backend_strides; 50822eb1385SJeremy L Thompson CeedInt layout[3] = {1, size, elem_size}; 50922eb1385SJeremy L Thompson 510dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetELayout(rstr, layout)); 51122eb1385SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_STRIDED) { 51222eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides)); 51322eb1385SJeremy L Thompson if (has_backend_strides) { 51422eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionSetLLayout(rstr, layout)); 51522eb1385SJeremy L Thompson } 51622eb1385SJeremy L Thompson } 51722eb1385SJeremy L Thompson } 518ff1e7120SSebastian Grimberg 519dce49693SSebastian Grimberg // Set up device offset/orientation arrays 520dce49693SSebastian Grimberg if (rstr_type != CEED_RESTRICTION_STRIDED) { 521ff1e7120SSebastian Grimberg switch (mem_type) { 522ff1e7120SSebastian Grimberg case CEED_MEM_HOST: { 523ff1e7120SSebastian Grimberg switch (copy_mode) { 524*a267acd1SJeremy L Thompson case CEED_COPY_VALUES: 525*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_offsets_owned)); 526*a267acd1SJeremy L Thompson memcpy(impl->h_offsets_owned, offsets, size * sizeof(CeedInt)); 527*a267acd1SJeremy L Thompson impl->h_offsets_borrowed = NULL; 528*a267acd1SJeremy L Thompson impl->h_offsets = impl->h_offsets_owned; 529*a267acd1SJeremy L Thompson break; 530ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 531*a267acd1SJeremy L Thompson impl->h_offsets_owned = (CeedInt *)offsets; 532*a267acd1SJeremy L Thompson impl->h_offsets_borrowed = NULL; 533*a267acd1SJeremy L Thompson impl->h_offsets = impl->h_offsets_owned; 534ff1e7120SSebastian Grimberg break; 535ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 536*a267acd1SJeremy L Thompson impl->h_offsets_owned = NULL; 537*a267acd1SJeremy L Thompson impl->h_offsets_borrowed = (CeedInt *)offsets; 538*a267acd1SJeremy L Thompson impl->h_offsets = impl->h_offsets_borrowed; 539ff1e7120SSebastian Grimberg break; 540ff1e7120SSebastian Grimberg } 541*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_owned, size * sizeof(CeedInt))); 542*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_offsets_owned, impl->h_offsets, size * sizeof(CeedInt), cudaMemcpyHostToDevice)); 543*a267acd1SJeremy L Thompson impl->d_offsets = impl->d_offsets_owned; 544*a267acd1SJeremy L Thompson if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, offsets)); 545dce49693SSebastian Grimberg } break; 546ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: { 547ff1e7120SSebastian Grimberg switch (copy_mode) { 548ff1e7120SSebastian Grimberg case CEED_COPY_VALUES: 549*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_owned, size * sizeof(CeedInt))); 550*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_offsets_owned, offsets, size * sizeof(CeedInt), cudaMemcpyDeviceToDevice)); 551*a267acd1SJeremy L Thompson impl->d_offsets_borrowed = NULL; 552*a267acd1SJeremy L Thompson impl->d_offsets = impl->d_offsets_owned; 553ff1e7120SSebastian Grimberg break; 554ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 555*a267acd1SJeremy L Thompson impl->d_offsets_owned = (CeedInt *)offsets; 556*a267acd1SJeremy L Thompson impl->d_offsets_borrowed = NULL; 557*a267acd1SJeremy L Thompson impl->d_offsets = impl->d_offsets_owned; 558ff1e7120SSebastian Grimberg break; 559ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 560*a267acd1SJeremy L Thompson impl->d_offsets_owned = NULL; 561*a267acd1SJeremy L Thompson impl->d_offsets_borrowed = (CeedInt *)offsets; 562*a267acd1SJeremy L Thompson impl->d_offsets = impl->d_offsets_borrowed; 563ff1e7120SSebastian Grimberg break; 564ff1e7120SSebastian Grimberg } 565*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_offsets_owned)); 566*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->h_offsets_owned, impl->d_offsets, size * sizeof(CeedInt), cudaMemcpyDeviceToHost)); 567*a267acd1SJeremy L Thompson impl->h_offsets = impl->h_offsets_owned; 568*a267acd1SJeremy L Thompson if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, offsets)); 569dce49693SSebastian Grimberg } break; 570ff1e7120SSebastian Grimberg } 571ff1e7120SSebastian Grimberg 572dce49693SSebastian Grimberg // Orientation data 573dce49693SSebastian Grimberg if (rstr_type == CEED_RESTRICTION_ORIENTED) { 574dce49693SSebastian Grimberg switch (mem_type) { 575dce49693SSebastian Grimberg case CEED_MEM_HOST: { 576dce49693SSebastian Grimberg switch (copy_mode) { 577*a267acd1SJeremy L Thompson case CEED_COPY_VALUES: 578*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_orients_owned)); 579*a267acd1SJeremy L Thompson memcpy(impl->h_orients_owned, orients, size * sizeof(bool)); 580*a267acd1SJeremy L Thompson impl->h_orients_borrowed = NULL; 581*a267acd1SJeremy L Thompson impl->h_orients = impl->h_orients_owned; 582*a267acd1SJeremy L Thompson break; 583dce49693SSebastian Grimberg case CEED_OWN_POINTER: 584*a267acd1SJeremy L Thompson impl->h_orients_owned = (bool *)orients; 585*a267acd1SJeremy L Thompson impl->h_orients_borrowed = NULL; 586*a267acd1SJeremy L Thompson impl->h_orients = impl->h_orients_owned; 587dce49693SSebastian Grimberg break; 588dce49693SSebastian Grimberg case CEED_USE_POINTER: 589*a267acd1SJeremy L Thompson impl->h_orients_owned = NULL; 590*a267acd1SJeremy L Thompson impl->h_orients_borrowed = (bool *)orients; 591*a267acd1SJeremy L Thompson impl->h_orients = impl->h_orients_borrowed; 592dce49693SSebastian Grimberg break; 593dce49693SSebastian Grimberg } 594*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_orients_owned, size * sizeof(bool))); 595*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_orients_owned, impl->h_orients, size * sizeof(bool), cudaMemcpyHostToDevice)); 596*a267acd1SJeremy L Thompson impl->d_orients = impl->d_orients_owned; 597dce49693SSebastian Grimberg } break; 598dce49693SSebastian Grimberg case CEED_MEM_DEVICE: { 599dce49693SSebastian Grimberg switch (copy_mode) { 600dce49693SSebastian Grimberg case CEED_COPY_VALUES: 601*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_orients_owned, size * sizeof(bool))); 602*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_orients_owned, orients, size * sizeof(bool), cudaMemcpyDeviceToDevice)); 603*a267acd1SJeremy L Thompson impl->d_orients_borrowed = NULL; 604*a267acd1SJeremy L Thompson impl->d_orients = impl->d_orients_owned; 605dce49693SSebastian Grimberg break; 606dce49693SSebastian Grimberg case CEED_OWN_POINTER: 607*a267acd1SJeremy L Thompson impl->d_orients_owned = (bool *)orients; 608*a267acd1SJeremy L Thompson impl->d_orients_borrowed = NULL; 609*a267acd1SJeremy L Thompson impl->d_orients = impl->d_orients_owned; 610dce49693SSebastian Grimberg break; 611dce49693SSebastian Grimberg case CEED_USE_POINTER: 612*a267acd1SJeremy L Thompson impl->d_orients_owned = NULL; 613*a267acd1SJeremy L Thompson impl->d_orients_borrowed = (bool *)orients; 614*a267acd1SJeremy L Thompson impl->d_orients = impl->d_orients_borrowed; 615dce49693SSebastian Grimberg break; 616dce49693SSebastian Grimberg } 617*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_orients_owned)); 618*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->h_orients_owned, impl->d_orients, size * sizeof(bool), cudaMemcpyDeviceToHost)); 619*a267acd1SJeremy L Thompson impl->h_orients = impl->h_orients_owned; 620dce49693SSebastian Grimberg } break; 621dce49693SSebastian Grimberg } 622dce49693SSebastian Grimberg } else if (rstr_type == CEED_RESTRICTION_CURL_ORIENTED) { 623dce49693SSebastian Grimberg switch (mem_type) { 624dce49693SSebastian Grimberg case CEED_MEM_HOST: { 625dce49693SSebastian Grimberg switch (copy_mode) { 626*a267acd1SJeremy L Thompson case CEED_COPY_VALUES: 627*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_owned)); 628*a267acd1SJeremy L Thompson memcpy(impl->h_curl_orients_owned, curl_orients, 3 * size * sizeof(CeedInt8)); 629*a267acd1SJeremy L Thompson impl->h_curl_orients_borrowed = NULL; 630*a267acd1SJeremy L Thompson impl->h_curl_orients = impl->h_curl_orients_owned; 631*a267acd1SJeremy L Thompson break; 632dce49693SSebastian Grimberg case CEED_OWN_POINTER: 633*a267acd1SJeremy L Thompson impl->h_curl_orients_owned = (CeedInt8 *)curl_orients; 634*a267acd1SJeremy L Thompson impl->h_curl_orients_borrowed = NULL; 635*a267acd1SJeremy L Thompson impl->h_curl_orients = impl->h_curl_orients_owned; 636dce49693SSebastian Grimberg break; 637dce49693SSebastian Grimberg case CEED_USE_POINTER: 638*a267acd1SJeremy L Thompson impl->h_curl_orients_owned = NULL; 639*a267acd1SJeremy L Thompson impl->h_curl_orients_borrowed = (CeedInt8 *)curl_orients; 640*a267acd1SJeremy L Thompson impl->h_curl_orients = impl->h_curl_orients_borrowed; 641dce49693SSebastian Grimberg break; 642dce49693SSebastian Grimberg } 643*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_curl_orients_owned, 3 * size * sizeof(CeedInt8))); 644*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_curl_orients_owned, impl->h_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyHostToDevice)); 645*a267acd1SJeremy L Thompson impl->d_curl_orients = impl->d_curl_orients_owned; 646dce49693SSebastian Grimberg } break; 647dce49693SSebastian Grimberg case CEED_MEM_DEVICE: { 648dce49693SSebastian Grimberg switch (copy_mode) { 649dce49693SSebastian Grimberg case CEED_COPY_VALUES: 650*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_curl_orients_owned, 3 * size * sizeof(CeedInt8))); 651*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->d_curl_orients_owned, curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyDeviceToDevice)); 652*a267acd1SJeremy L Thompson impl->d_curl_orients_borrowed = NULL; 653*a267acd1SJeremy L Thompson impl->d_curl_orients = impl->d_curl_orients_owned; 654dce49693SSebastian Grimberg break; 655dce49693SSebastian Grimberg case CEED_OWN_POINTER: 656*a267acd1SJeremy L Thompson impl->d_curl_orients_owned = (CeedInt8 *)curl_orients; 657*a267acd1SJeremy L Thompson impl->d_curl_orients_borrowed = NULL; 658*a267acd1SJeremy L Thompson impl->d_curl_orients = impl->d_curl_orients_owned; 659dce49693SSebastian Grimberg break; 660dce49693SSebastian Grimberg case CEED_USE_POINTER: 661*a267acd1SJeremy L Thompson impl->d_curl_orients_owned = NULL; 662*a267acd1SJeremy L Thompson impl->d_curl_orients_borrowed = (CeedInt8 *)curl_orients; 663*a267acd1SJeremy L Thompson impl->d_curl_orients = impl->d_curl_orients_borrowed; 664dce49693SSebastian Grimberg break; 665dce49693SSebastian Grimberg } 666*a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_owned)); 667*a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->h_curl_orients_owned, impl->d_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyDeviceToHost)); 668*a267acd1SJeremy L Thompson impl->h_curl_orients = impl->h_curl_orients_owned; 669dce49693SSebastian Grimberg } break; 670dce49693SSebastian Grimberg } 671dce49693SSebastian Grimberg } 672dce49693SSebastian Grimberg } 673ca735530SJeremy L Thompson 674ff1e7120SSebastian Grimberg // Register backend functions 675dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Apply", CeedElemRestrictionApply_Cuda)); 676dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnsigned", CeedElemRestrictionApplyUnsigned_Cuda)); 677dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnoriented", CeedElemRestrictionApplyUnoriented_Cuda)); 678dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOffsets", CeedElemRestrictionGetOffsets_Cuda)); 679dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOrientations", CeedElemRestrictionGetOrientations_Cuda)); 680dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetCurlOrientations", CeedElemRestrictionGetCurlOrientations_Cuda)); 681dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Destroy", CeedElemRestrictionDestroy_Cuda)); 682ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 683ff1e7120SSebastian Grimberg } 684ff1e7120SSebastian Grimberg 685ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 686