15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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)); 35*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 36cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 37cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 38cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride)); 39*b20a4af9SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_POINTS) { 40*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr, &elem_size)); 41*b20a4af9SJeremy L Thompson } else { 42*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 43*b20a4af9SJeremy L Thompson } 44cf8cbdd6SSebastian Grimberg is_deterministic = impl->d_l_vec_indices != NULL; 45cf8cbdd6SSebastian Grimberg 46cf8cbdd6SSebastian Grimberg // Compile CUDA kernels 47cf8cbdd6SSebastian Grimberg switch (rstr_type) { 48cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 49cf8cbdd6SSebastian Grimberg bool has_backend_strides; 50509d4af6SJeremy L Thompson CeedInt strides[3] = {1, num_elem * elem_size, elem_size}; 51cf8cbdd6SSebastian Grimberg 52cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides)); 53cf8cbdd6SSebastian Grimberg if (!has_backend_strides) { 5456c48462SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetStrides(rstr, strides)); 55cf8cbdd6SSebastian Grimberg } 56cf8cbdd6SSebastian Grimberg 57cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-strided.h", &restriction_kernel_path)); 58cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 59cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 60cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 61cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 62cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_STRIDE_NODES", strides[0], "RSTR_STRIDE_COMP", strides[1], "RSTR_STRIDE_ELEM", 63cf8cbdd6SSebastian Grimberg strides[2])); 64cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedNoTranspose", &impl->ApplyNoTranspose)); 65cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "StridedTranspose", &impl->ApplyTranspose)); 66cf8cbdd6SSebastian Grimberg } break; 67*b20a4af9SJeremy L Thompson case CEED_RESTRICTION_POINTS: 68cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 69cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &restriction_kernel_path)); 70cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 71cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedLoadSourceToBuffer(ceed, restriction_kernel_path, &restriction_kernel_source)); 72cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 73cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 74cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 75cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 76cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyNoTranspose)); 77cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyTranspose)); 78cf8cbdd6SSebastian Grimberg } break; 79cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 8022070f95SJeremy L Thompson const char *offset_kernel_path; 81509d4af6SJeremy L Thompson char **file_paths = NULL; 82509d4af6SJeremy L Thompson CeedInt num_file_paths = 0; 83cf8cbdd6SSebastian Grimberg 84cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-oriented.h", &restriction_kernel_path)); 85cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 86509d4af6SJeremy L Thompson CeedCallBackend(CeedLoadSourceAndInitializeBuffer(ceed, restriction_kernel_path, &num_file_paths, &file_paths, &restriction_kernel_source)); 87cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &offset_kernel_path)); 88509d4af6SJeremy L Thompson CeedCallBackend(CeedLoadSourceToInitializedBuffer(ceed, offset_kernel_path, &num_file_paths, &file_paths, &restriction_kernel_source)); 89cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 90cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 91cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 92cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 93cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedNoTranspose", &impl->ApplyNoTranspose)); 94cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnsignedNoTranspose)); 95cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OrientedTranspose", &impl->ApplyTranspose)); 96cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnsignedTranspose)); 97509d4af6SJeremy L Thompson // Cleanup 98cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&offset_kernel_path)); 99509d4af6SJeremy L Thompson for (CeedInt i = 0; i < num_file_paths; i++) CeedCall(CeedFree(&file_paths[i])); 100509d4af6SJeremy L Thompson CeedCall(CeedFree(&file_paths)); 101cf8cbdd6SSebastian Grimberg } break; 102cf8cbdd6SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 10322070f95SJeremy L Thompson const char *offset_kernel_path; 104509d4af6SJeremy L Thompson char **file_paths = NULL; 105509d4af6SJeremy L Thompson CeedInt num_file_paths = 0; 106cf8cbdd6SSebastian Grimberg 107cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-curl-oriented.h", &restriction_kernel_path)); 108cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source -----\n"); 109509d4af6SJeremy L Thompson CeedCallBackend(CeedLoadSourceAndInitializeBuffer(ceed, restriction_kernel_path, &num_file_paths, &file_paths, &restriction_kernel_source)); 110cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-restriction-offset.h", &offset_kernel_path)); 111509d4af6SJeremy L Thompson CeedCallBackend(CeedLoadSourceToInitializedBuffer(ceed, offset_kernel_path, &num_file_paths, &file_paths, &restriction_kernel_source)); 112cf8cbdd6SSebastian Grimberg CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Restriction Kernel Source Complete! -----\n"); 113cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedCompile_Cuda(ceed, restriction_kernel_source, &impl->module, 6, "RSTR_ELEM_SIZE", elem_size, "RSTR_NUM_ELEM", num_elem, 114cf8cbdd6SSebastian Grimberg "RSTR_NUM_COMP", num_comp, "RSTR_NUM_NODES", impl->num_nodes, "RSTR_COMP_STRIDE", comp_stride, 115cf8cbdd6SSebastian Grimberg "USE_DETERMINISTIC", is_deterministic ? 1 : 0)); 116cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedNoTranspose", &impl->ApplyNoTranspose)); 117cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedNoTranspose", &impl->ApplyUnsignedNoTranspose)); 118cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetNoTranspose", &impl->ApplyUnorientedNoTranspose)); 119cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedTranspose", &impl->ApplyTranspose)); 120cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "CurlOrientedUnsignedTranspose", &impl->ApplyUnsignedTranspose)); 121cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedGetKernel_Cuda(ceed, impl->module, "OffsetTranspose", &impl->ApplyUnorientedTranspose)); 122509d4af6SJeremy L Thompson // Cleanup 123cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&offset_kernel_path)); 124509d4af6SJeremy L Thompson for (CeedInt i = 0; i < num_file_paths; i++) CeedCall(CeedFree(&file_paths[i])); 125509d4af6SJeremy L Thompson CeedCall(CeedFree(&file_paths)); 126cf8cbdd6SSebastian Grimberg } break; 127cf8cbdd6SSebastian Grimberg } 128cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_path)); 129cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedFree(&restriction_kernel_source)); 130cf8cbdd6SSebastian Grimberg return CEED_ERROR_SUCCESS; 131cf8cbdd6SSebastian Grimberg } 132cf8cbdd6SSebastian Grimberg 133cf8cbdd6SSebastian Grimberg //------------------------------------------------------------------------------ 134dce49693SSebastian Grimberg // Core apply restriction code 135ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 136dce49693SSebastian Grimberg static inline int CeedElemRestrictionApply_Cuda_Core(CeedElemRestriction rstr, CeedTransposeMode t_mode, bool use_signs, bool use_orients, 137dce49693SSebastian Grimberg CeedVector u, CeedVector v, CeedRequest *request) { 138ff1e7120SSebastian Grimberg Ceed ceed; 139dce49693SSebastian Grimberg CeedRestrictionType rstr_type; 140ff1e7120SSebastian Grimberg const CeedScalar *d_u; 141ff1e7120SSebastian Grimberg CeedScalar *d_v; 142ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 143ca735530SJeremy L Thompson 144dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 145dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 146dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 147cf8cbdd6SSebastian Grimberg 148cf8cbdd6SSebastian Grimberg // Assemble kernel if needed 149cf8cbdd6SSebastian Grimberg if (!impl->module) { 150cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetupCompile_Cuda(rstr)); 151cf8cbdd6SSebastian Grimberg } 152ca735530SJeremy L Thompson 153ca735530SJeremy L Thompson // Get vectors 154ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 155ff1e7120SSebastian Grimberg if (t_mode == CEED_TRANSPOSE) { 156ff1e7120SSebastian Grimberg // Sum into for transpose mode, e-vec to l-vec 157ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArray(v, CEED_MEM_DEVICE, &d_v)); 158ff1e7120SSebastian Grimberg } else { 159ff1e7120SSebastian Grimberg // Overwrite for notranspose mode, l-vec to e-vec 160ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 161ff1e7120SSebastian Grimberg } 162ff1e7120SSebastian Grimberg 163ff1e7120SSebastian Grimberg // Restrict 164ff1e7120SSebastian Grimberg if (t_mode == CEED_NOTRANSPOSE) { 165ff1e7120SSebastian Grimberg // L-vector -> E-vector 166cf8cbdd6SSebastian Grimberg CeedInt elem_size; 167cf8cbdd6SSebastian Grimberg 168cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 169dce49693SSebastian Grimberg const CeedInt block_size = elem_size < 1024 ? (elem_size > 32 ? elem_size : 32) : 1024; 170cf8cbdd6SSebastian Grimberg const CeedInt grid = CeedDivUpInt(impl->num_nodes, block_size); 17158549094SSebastian Grimberg 172dce49693SSebastian Grimberg switch (rstr_type) { 173dce49693SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 174cf8cbdd6SSebastian Grimberg void *args[] = {&d_u, &d_v}; 17558549094SSebastian Grimberg 176cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 177dce49693SSebastian Grimberg } break; 178*b20a4af9SJeremy L Thompson case CEED_RESTRICTION_POINTS: 179dce49693SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 180a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 181dce49693SSebastian Grimberg 182cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 183dce49693SSebastian Grimberg } break; 184dce49693SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 185dce49693SSebastian Grimberg if (use_signs) { 186a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v}; 187dce49693SSebastian Grimberg 188cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 189dce49693SSebastian Grimberg } else { 190a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 191dce49693SSebastian Grimberg 192cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args)); 193dce49693SSebastian Grimberg } 194dce49693SSebastian Grimberg } break; 195dce49693SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 196dce49693SSebastian Grimberg if (use_signs && use_orients) { 197a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 198dce49693SSebastian Grimberg 199cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyNoTranspose, grid, block_size, args)); 200dce49693SSebastian Grimberg } else if (use_orients) { 201a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 202dce49693SSebastian Grimberg 203cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedNoTranspose, grid, block_size, args)); 204dce49693SSebastian Grimberg } else { 205a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 206dce49693SSebastian Grimberg 207cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedNoTranspose, grid, block_size, args)); 208dce49693SSebastian Grimberg } 209dce49693SSebastian Grimberg } break; 210ff1e7120SSebastian Grimberg } 211ff1e7120SSebastian Grimberg } else { 212ff1e7120SSebastian Grimberg // E-vector -> L-vector 213cf8cbdd6SSebastian Grimberg const bool is_deterministic = impl->d_l_vec_indices != NULL; 214dce49693SSebastian Grimberg const CeedInt block_size = 32; 215cf8cbdd6SSebastian Grimberg const CeedInt grid = CeedDivUpInt(impl->num_nodes, block_size); 216ca735530SJeremy L Thompson 217dce49693SSebastian Grimberg switch (rstr_type) { 218dce49693SSebastian Grimberg case CEED_RESTRICTION_STRIDED: { 219cf8cbdd6SSebastian Grimberg void *args[] = {&d_u, &d_v}; 220dce49693SSebastian Grimberg 221cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 222dce49693SSebastian Grimberg } break; 223*b20a4af9SJeremy L Thompson case CEED_RESTRICTION_POINTS: 224dce49693SSebastian Grimberg case CEED_RESTRICTION_STANDARD: { 225cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 226a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 22758549094SSebastian Grimberg 228cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 229ff1e7120SSebastian Grimberg } else { 23058549094SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 23158549094SSebastian Grimberg 232cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 23358549094SSebastian Grimberg } 234dce49693SSebastian Grimberg } break; 235dce49693SSebastian Grimberg case CEED_RESTRICTION_ORIENTED: { 236dce49693SSebastian Grimberg if (use_signs) { 237cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 238a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_orients, &d_u, &d_v}; 23958549094SSebastian Grimberg 240cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 241dce49693SSebastian Grimberg } else { 2427aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_orients, &d_u, &d_v}; 2437aa91133SSebastian Grimberg 244cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2457aa91133SSebastian Grimberg } 2467aa91133SSebastian Grimberg } else { 247cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 248a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 249dce49693SSebastian Grimberg 250cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 251dce49693SSebastian Grimberg } else { 252dce49693SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 253dce49693SSebastian Grimberg 254cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 255dce49693SSebastian Grimberg } 256dce49693SSebastian Grimberg } 257dce49693SSebastian Grimberg } break; 258dce49693SSebastian Grimberg case CEED_RESTRICTION_CURL_ORIENTED: { 259dce49693SSebastian Grimberg if (use_signs && use_orients) { 260cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 261a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 262dce49693SSebastian Grimberg 263cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2647aa91133SSebastian Grimberg } else { 2657aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v}; 2667aa91133SSebastian Grimberg 267cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyTranspose, grid, block_size, args)); 2687aa91133SSebastian Grimberg } 269dce49693SSebastian Grimberg } else if (use_orients) { 270cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 271a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &impl->d_curl_orients, &d_u, &d_v}; 272dce49693SSebastian Grimberg 273cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 274dce49693SSebastian Grimberg } else { 2757aa91133SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &impl->d_curl_orients, &d_u, &d_v}; 2767aa91133SSebastian Grimberg 277cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnsignedTranspose, grid, block_size, args)); 2787aa91133SSebastian Grimberg } 2797aa91133SSebastian Grimberg } else { 280cf8cbdd6SSebastian Grimberg if (!is_deterministic) { 281a267acd1SJeremy L Thompson void *args[] = {&impl->d_offsets, &d_u, &d_v}; 282dce49693SSebastian Grimberg 283cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args)); 284dce49693SSebastian Grimberg } else { 285dce49693SSebastian Grimberg void *args[] = {&impl->d_l_vec_indices, &impl->d_t_indices, &impl->d_t_offsets, &d_u, &d_v}; 286dce49693SSebastian Grimberg 287cf8cbdd6SSebastian Grimberg CeedCallBackend(CeedRunKernel_Cuda(ceed, impl->ApplyUnorientedTranspose, grid, block_size, args)); 288dce49693SSebastian Grimberg } 289dce49693SSebastian Grimberg } 290dce49693SSebastian Grimberg } break; 291ff1e7120SSebastian Grimberg } 292ff1e7120SSebastian Grimberg } 293ff1e7120SSebastian Grimberg 294ff1e7120SSebastian Grimberg if (request != CEED_REQUEST_IMMEDIATE && request != CEED_REQUEST_ORDERED) *request = NULL; 295ff1e7120SSebastian Grimberg 296ff1e7120SSebastian Grimberg // Restore arrays 297ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 298ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 299ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 300ff1e7120SSebastian Grimberg } 301ff1e7120SSebastian Grimberg 302ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 303dce49693SSebastian Grimberg // Apply restriction 304dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 305dce49693SSebastian Grimberg static int CeedElemRestrictionApply_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, CeedRequest *request) { 306dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, true, true, u, v, request); 307dce49693SSebastian Grimberg } 308dce49693SSebastian Grimberg 309dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 310dce49693SSebastian Grimberg // Apply unsigned restriction 311dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 312dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnsigned_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, 313dce49693SSebastian Grimberg CeedRequest *request) { 314dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, true, u, v, request); 315dce49693SSebastian Grimberg } 316dce49693SSebastian Grimberg 317dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 318dce49693SSebastian Grimberg // Apply unoriented restriction 319dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 320dce49693SSebastian Grimberg static int CeedElemRestrictionApplyUnoriented_Cuda(CeedElemRestriction rstr, CeedTransposeMode t_mode, CeedVector u, CeedVector v, 321dce49693SSebastian Grimberg CeedRequest *request) { 322dce49693SSebastian Grimberg return CeedElemRestrictionApply_Cuda_Core(rstr, t_mode, false, false, u, v, request); 323dce49693SSebastian Grimberg } 324dce49693SSebastian Grimberg 325dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 326ff1e7120SSebastian Grimberg // Get offsets 327ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 328ff1e7120SSebastian Grimberg static int CeedElemRestrictionGetOffsets_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt **offsets) { 329ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 330*b20a4af9SJeremy L Thompson CeedRestrictionType rstr_type; 331ff1e7120SSebastian Grimberg 332ca735530SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 333*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 334ff1e7120SSebastian Grimberg switch (mem_type) { 335ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 336*b20a4af9SJeremy L Thompson *offsets = rstr_type == CEED_RESTRICTION_POINTS ? impl->h_offsets_at_points : impl->h_offsets; 337ff1e7120SSebastian Grimberg break; 338ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 339*b20a4af9SJeremy L Thompson *offsets = rstr_type == CEED_RESTRICTION_POINTS ? impl->d_offsets_at_points : impl->d_offsets; 340ff1e7120SSebastian Grimberg break; 341ff1e7120SSebastian Grimberg } 342ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 343ff1e7120SSebastian Grimberg } 344ff1e7120SSebastian Grimberg 345ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 346dce49693SSebastian Grimberg // Get orientations 347dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 348dce49693SSebastian Grimberg static int CeedElemRestrictionGetOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const bool **orients) { 349dce49693SSebastian Grimberg CeedElemRestriction_Cuda *impl; 350dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 351dce49693SSebastian Grimberg 352dce49693SSebastian Grimberg switch (mem_type) { 353dce49693SSebastian Grimberg case CEED_MEM_HOST: 354dce49693SSebastian Grimberg *orients = impl->h_orients; 355dce49693SSebastian Grimberg break; 356dce49693SSebastian Grimberg case CEED_MEM_DEVICE: 357dce49693SSebastian Grimberg *orients = impl->d_orients; 358dce49693SSebastian Grimberg break; 359dce49693SSebastian Grimberg } 360dce49693SSebastian Grimberg return CEED_ERROR_SUCCESS; 361dce49693SSebastian Grimberg } 362dce49693SSebastian Grimberg 363dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 364dce49693SSebastian Grimberg // Get curl-conforming orientations 365dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 366dce49693SSebastian Grimberg static int CeedElemRestrictionGetCurlOrientations_Cuda(CeedElemRestriction rstr, CeedMemType mem_type, const CeedInt8 **curl_orients) { 367dce49693SSebastian Grimberg CeedElemRestriction_Cuda *impl; 368dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 369dce49693SSebastian Grimberg 370dce49693SSebastian Grimberg switch (mem_type) { 371dce49693SSebastian Grimberg case CEED_MEM_HOST: 372dce49693SSebastian Grimberg *curl_orients = impl->h_curl_orients; 373dce49693SSebastian Grimberg break; 374dce49693SSebastian Grimberg case CEED_MEM_DEVICE: 375dce49693SSebastian Grimberg *curl_orients = impl->d_curl_orients; 376dce49693SSebastian Grimberg break; 377dce49693SSebastian Grimberg } 378dce49693SSebastian Grimberg return CEED_ERROR_SUCCESS; 379dce49693SSebastian Grimberg } 380dce49693SSebastian Grimberg 381dce49693SSebastian Grimberg //------------------------------------------------------------------------------ 382*b20a4af9SJeremy L Thompson // Get offset for padded AtPoints E-layout 383*b20a4af9SJeremy L Thompson //------------------------------------------------------------------------------ 384*b20a4af9SJeremy L Thompson static int CeedElemRestrictionGetAtPointsElementOffset_Cuda(CeedElemRestriction rstr, CeedInt elem, CeedSize *elem_offset) { 385*b20a4af9SJeremy L Thompson CeedInt layout[3]; 386*b20a4af9SJeremy L Thompson 387*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetELayout(rstr, layout)); 388*b20a4af9SJeremy L Thompson *elem_offset = 0 * layout[0] + 0 * layout[1] + elem * layout[2]; 389*b20a4af9SJeremy L Thompson return CEED_ERROR_SUCCESS; 390*b20a4af9SJeremy L Thompson } 391*b20a4af9SJeremy L Thompson 392*b20a4af9SJeremy L Thompson //------------------------------------------------------------------------------ 393ff1e7120SSebastian Grimberg // Destroy restriction 394ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 395dce49693SSebastian Grimberg static int CeedElemRestrictionDestroy_Cuda(CeedElemRestriction rstr) { 396ff1e7120SSebastian Grimberg Ceed ceed; 397ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 398ca735530SJeremy L Thompson 399dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 400dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 401cf8cbdd6SSebastian Grimberg if (impl->module) { 402ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cuModuleUnload(impl->module)); 403cf8cbdd6SSebastian Grimberg } 404a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_offsets_owned)); 405f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_offsets_owned)); 406081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_t_offsets)); 407081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_t_indices)); 408081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_l_vec_indices)); 409a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_orients_owned)); 410f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaFree((bool *)impl->d_orients_owned)); 411a267acd1SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_curl_orients_owned)); 412f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt8 *)impl->d_curl_orients_owned)); 413*b20a4af9SJeremy L Thompson CeedCallBackend(CeedFree(&impl->h_offsets_at_points_owned)); 414*b20a4af9SJeremy L Thompson CeedCallCuda(ceed, cudaFree((CeedInt *)impl->d_offsets_at_points_owned)); 415ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 416ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 417ff1e7120SSebastian Grimberg } 418ff1e7120SSebastian Grimberg 419ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 420ff1e7120SSebastian Grimberg // Create transpose offsets and indices 421ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 422*b20a4af9SJeremy L Thompson static int CeedElemRestrictionOffset_Cuda(const CeedElemRestriction rstr, const CeedInt elem_size, const CeedInt *indices) { 423ff1e7120SSebastian Grimberg Ceed ceed; 424ca735530SJeremy L Thompson bool *is_node; 425ff1e7120SSebastian Grimberg CeedSize l_size; 426*b20a4af9SJeremy L Thompson CeedInt num_elem, num_comp, num_nodes = 0; 427ca735530SJeremy L Thompson CeedInt *ind_to_offset, *l_vec_indices, *t_offsets, *t_indices; 428*b20a4af9SJeremy L Thompson CeedRestrictionType rstr_type; 429ca735530SJeremy L Thompson CeedElemRestriction_Cuda *impl; 430ca735530SJeremy L Thompson 431dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 432dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetData(rstr, &impl)); 433dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 434*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 435dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size)); 436dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 437ca735530SJeremy L Thompson const CeedInt size_indices = num_elem * elem_size; 438ff1e7120SSebastian Grimberg 439ff1e7120SSebastian Grimberg // Count num_nodes 440ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &is_node)); 441ca735530SJeremy L Thompson 442ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < size_indices; i++) is_node[indices[i]] = 1; 443ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < l_size; i++) num_nodes += is_node[i]; 444ff1e7120SSebastian Grimberg impl->num_nodes = num_nodes; 445ff1e7120SSebastian Grimberg 446ff1e7120SSebastian Grimberg // L-vector offsets array 447ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(l_size, &ind_to_offset)); 448ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(num_nodes, &l_vec_indices)); 449ca735530SJeremy L Thompson for (CeedInt i = 0, j = 0; i < l_size; i++) { 450ff1e7120SSebastian Grimberg if (is_node[i]) { 451ff1e7120SSebastian Grimberg l_vec_indices[j] = i; 452ff1e7120SSebastian Grimberg ind_to_offset[i] = j++; 453ff1e7120SSebastian Grimberg } 454ff1e7120SSebastian Grimberg } 455ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&is_node)); 456ff1e7120SSebastian Grimberg 457ff1e7120SSebastian Grimberg // Compute transpose offsets and indices 458ff1e7120SSebastian Grimberg const CeedInt size_offsets = num_nodes + 1; 459ca735530SJeremy L Thompson 460ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(size_offsets, &t_offsets)); 461ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(size_indices, &t_indices)); 462ff1e7120SSebastian Grimberg // Count node multiplicity 463ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 464ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) ++t_offsets[ind_to_offset[indices[elem_size * e + i]] + 1]; 465ff1e7120SSebastian Grimberg } 466ff1e7120SSebastian Grimberg // Convert to running sum 467ff1e7120SSebastian Grimberg for (CeedInt i = 1; i < size_offsets; ++i) t_offsets[i] += t_offsets[i - 1]; 468ff1e7120SSebastian Grimberg // List all E-vec indices associated with L-vec node 469ff1e7120SSebastian Grimberg for (CeedInt e = 0; e < num_elem; ++e) { 470ff1e7120SSebastian Grimberg for (CeedInt i = 0; i < elem_size; ++i) { 471ff1e7120SSebastian Grimberg const CeedInt lid = elem_size * e + i; 472ff1e7120SSebastian Grimberg const CeedInt gid = indices[lid]; 473ca735530SJeremy L Thompson 474ff1e7120SSebastian Grimberg t_indices[t_offsets[ind_to_offset[gid]]++] = lid; 475ff1e7120SSebastian Grimberg } 476ff1e7120SSebastian Grimberg } 477ff1e7120SSebastian Grimberg // Reset running sum 478ff1e7120SSebastian Grimberg for (int i = size_offsets - 1; i > 0; --i) t_offsets[i] = t_offsets[i - 1]; 479ff1e7120SSebastian Grimberg t_offsets[0] = 0; 480ff1e7120SSebastian Grimberg 481ff1e7120SSebastian Grimberg // Copy data to device 482ff1e7120SSebastian Grimberg // -- L-vector indices 483ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_l_vec_indices, num_nodes * sizeof(CeedInt))); 484081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_l_vec_indices, l_vec_indices, num_nodes * sizeof(CeedInt), cudaMemcpyHostToDevice)); 485ff1e7120SSebastian Grimberg // -- Transpose offsets 486ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_offsets, size_offsets * sizeof(CeedInt))); 487081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_t_offsets, t_offsets, size_offsets * sizeof(CeedInt), cudaMemcpyHostToDevice)); 488ff1e7120SSebastian Grimberg // -- Transpose indices 489ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_t_indices, size_indices * sizeof(CeedInt))); 490081aa29dSJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_t_indices, t_indices, size_indices * sizeof(CeedInt), cudaMemcpyHostToDevice)); 491ff1e7120SSebastian Grimberg 492ff1e7120SSebastian Grimberg // Cleanup 493ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&ind_to_offset)); 494ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&l_vec_indices)); 495ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_offsets)); 496ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&t_indices)); 497ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 498ff1e7120SSebastian Grimberg } 499ff1e7120SSebastian Grimberg 500ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 501ff1e7120SSebastian Grimberg // Create restriction 502ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 503a267acd1SJeremy L Thompson int CeedElemRestrictionCreate_Cuda(CeedMemType mem_type, CeedCopyMode copy_mode, const CeedInt *offsets, const bool *orients, 504dce49693SSebastian Grimberg const CeedInt8 *curl_orients, CeedElemRestriction rstr) { 505ca735530SJeremy L Thompson Ceed ceed, ceed_parent; 506dce49693SSebastian Grimberg bool is_deterministic; 507*b20a4af9SJeremy L Thompson CeedInt num_elem, num_comp, elem_size; 508ca735530SJeremy L Thompson CeedRestrictionType rstr_type; 509ff1e7120SSebastian Grimberg CeedElemRestriction_Cuda *impl; 510ca735530SJeremy L Thompson 511dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed)); 512ca735530SJeremy L Thompson CeedCallBackend(CeedGetParent(ceed, &ceed_parent)); 513ca735530SJeremy L Thompson CeedCallBackend(CeedIsDeterministic(ceed_parent, &is_deterministic)); 514dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 515*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 516dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 51722eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 518*b20a4af9SJeremy L Thompson // Use max number of points as elem size for AtPoints restrictions 519*b20a4af9SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_POINTS) { 520*b20a4af9SJeremy L Thompson CeedInt max_points = 0; 521*b20a4af9SJeremy L Thompson 522*b20a4af9SJeremy L Thompson for (CeedInt i = 0; i < num_elem; i++) { 523*b20a4af9SJeremy L Thompson max_points = CeedIntMax(max_points, offsets[i + 1] - offsets[i]); 524*b20a4af9SJeremy L Thompson } 525*b20a4af9SJeremy L Thompson elem_size = max_points; 526*b20a4af9SJeremy L Thompson } 527ca735530SJeremy L Thompson const CeedInt size = num_elem * elem_size; 528ff1e7120SSebastian Grimberg 529dce49693SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 530dce49693SSebastian Grimberg impl->num_nodes = size; 531dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetData(rstr, impl)); 53222eb1385SJeremy L Thompson 53322eb1385SJeremy L Thompson // Set layouts 53422eb1385SJeremy L Thompson { 53522eb1385SJeremy L Thompson bool has_backend_strides; 53622eb1385SJeremy L Thompson CeedInt layout[3] = {1, size, elem_size}; 53722eb1385SJeremy L Thompson 538dce49693SSebastian Grimberg CeedCallBackend(CeedElemRestrictionSetELayout(rstr, layout)); 53922eb1385SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_STRIDED) { 54022eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionHasBackendStrides(rstr, &has_backend_strides)); 54122eb1385SJeremy L Thompson if (has_backend_strides) { 54222eb1385SJeremy L Thompson CeedCallBackend(CeedElemRestrictionSetLLayout(rstr, layout)); 54322eb1385SJeremy L Thompson } 54422eb1385SJeremy L Thompson } 54522eb1385SJeremy L Thompson } 546ff1e7120SSebastian Grimberg 547*b20a4af9SJeremy L Thompson // Pad AtPoints indices 548*b20a4af9SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_POINTS) { 549*b20a4af9SJeremy L Thompson CeedSize offsets_len = elem_size * num_elem, at_points_size = num_elem + 1; 550*b20a4af9SJeremy L Thompson CeedInt max_points = elem_size, *offsets_padded; 551*b20a4af9SJeremy L Thompson 552*b20a4af9SJeremy L Thompson CeedCheck(mem_type == CEED_MEM_HOST, ceed, CEED_ERROR_BACKEND, "only MemType Host supported when creating AtPoints restriction"); 553*b20a4af9SJeremy L Thompson CeedCallBackend(CeedMalloc(offsets_len, &offsets_padded)); 554*b20a4af9SJeremy L Thompson for (CeedInt i = 0; i < num_elem; i++) { 555*b20a4af9SJeremy L Thompson CeedInt num_points = offsets[i + 1] - offsets[i]; 556*b20a4af9SJeremy L Thompson 557*b20a4af9SJeremy L Thompson at_points_size += num_points; 558*b20a4af9SJeremy L Thompson // -- Copy all points in element 559*b20a4af9SJeremy L Thompson for (CeedInt j = 0; j < num_points; j++) { 560*b20a4af9SJeremy L Thompson offsets_padded[i * max_points + j] = offsets[offsets[i] + j]; 561*b20a4af9SJeremy L Thompson } 562*b20a4af9SJeremy L Thompson // -- Replicate out last point in element 563*b20a4af9SJeremy L Thompson for (CeedInt j = num_points; j < max_points; j++) { 564*b20a4af9SJeremy L Thompson offsets_padded[i * max_points + j] = offsets[offsets[i] + num_points - 1]; 565*b20a4af9SJeremy L Thompson } 566*b20a4af9SJeremy L Thompson } 567*b20a4af9SJeremy L Thompson CeedCallBackend(CeedSetHostCeedIntArray(offsets, copy_mode, at_points_size, &impl->h_offsets_at_points_owned, &impl->h_offsets_at_points_borrowed, 568*b20a4af9SJeremy L Thompson &impl->h_offsets_at_points)); 569*b20a4af9SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_at_points_owned, at_points_size * sizeof(CeedInt))); 570*b20a4af9SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt **)impl->d_offsets_at_points_owned, impl->h_offsets_at_points, at_points_size * sizeof(CeedInt), 571*b20a4af9SJeremy L Thompson cudaMemcpyHostToDevice)); 572*b20a4af9SJeremy L Thompson impl->d_offsets_at_points = (CeedInt *)impl->d_offsets_at_points_owned; 573*b20a4af9SJeremy L Thompson 574*b20a4af9SJeremy L Thompson // -- Use padded offsets for the rest of the setup 575*b20a4af9SJeremy L Thompson offsets = (const CeedInt *)offsets_padded; 576*b20a4af9SJeremy L Thompson copy_mode = CEED_OWN_POINTER; 577*b20a4af9SJeremy L Thompson CeedCallBackend(CeedElemRestrictionSetAtPointsEVectorSize(rstr, at_points_size * num_comp)); 578*b20a4af9SJeremy L Thompson } 579*b20a4af9SJeremy L Thompson 580dce49693SSebastian Grimberg // Set up device offset/orientation arrays 581dce49693SSebastian Grimberg if (rstr_type != CEED_RESTRICTION_STRIDED) { 582ff1e7120SSebastian Grimberg switch (mem_type) { 583ff1e7120SSebastian Grimberg case CEED_MEM_HOST: { 584f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetHostCeedIntArray(offsets, copy_mode, size, &impl->h_offsets_owned, &impl->h_offsets_borrowed, &impl->h_offsets)); 585a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_offsets_owned, size * sizeof(CeedInt))); 586f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->d_offsets_owned, impl->h_offsets, size * sizeof(CeedInt), cudaMemcpyHostToDevice)); 587f5d1e504SJeremy L Thompson impl->d_offsets = (CeedInt *)impl->d_offsets_owned; 588*b20a4af9SJeremy L Thompson if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, elem_size, offsets)); 589dce49693SSebastian Grimberg } break; 590ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: { 591f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetDeviceCeedIntArray_Cuda(ceed, offsets, copy_mode, size, &impl->d_offsets_owned, &impl->d_offsets_borrowed, 592f5d1e504SJeremy L Thompson (const CeedInt **)&impl->d_offsets)); 593a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_offsets_owned)); 594f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((CeedInt *)impl->h_offsets_owned, impl->d_offsets, size * sizeof(CeedInt), cudaMemcpyDeviceToHost)); 595a267acd1SJeremy L Thompson impl->h_offsets = impl->h_offsets_owned; 596*b20a4af9SJeremy L Thompson if (is_deterministic) CeedCallBackend(CeedElemRestrictionOffset_Cuda(rstr, elem_size, offsets)); 597dce49693SSebastian Grimberg } break; 598ff1e7120SSebastian Grimberg } 599ff1e7120SSebastian Grimberg 600dce49693SSebastian Grimberg // Orientation data 601dce49693SSebastian Grimberg if (rstr_type == CEED_RESTRICTION_ORIENTED) { 602dce49693SSebastian Grimberg switch (mem_type) { 603dce49693SSebastian Grimberg case CEED_MEM_HOST: { 604f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetHostBoolArray(orients, copy_mode, size, &impl->h_orients_owned, &impl->h_orients_borrowed, &impl->h_orients)); 605a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_orients_owned, size * sizeof(bool))); 606f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((bool *)impl->d_orients_owned, impl->h_orients, size * sizeof(bool), cudaMemcpyHostToDevice)); 607a267acd1SJeremy L Thompson impl->d_orients = impl->d_orients_owned; 608dce49693SSebastian Grimberg } break; 609dce49693SSebastian Grimberg case CEED_MEM_DEVICE: { 610f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetDeviceBoolArray_Cuda(ceed, orients, copy_mode, size, &impl->d_orients_owned, &impl->d_orients_borrowed, 611f5d1e504SJeremy L Thompson (const bool **)&impl->d_orients)); 612a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(size, &impl->h_orients_owned)); 613f5d1e504SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy((bool *)impl->h_orients_owned, impl->d_orients, size * sizeof(bool), cudaMemcpyDeviceToHost)); 614a267acd1SJeremy L Thompson impl->h_orients = impl->h_orients_owned; 615dce49693SSebastian Grimberg } break; 616dce49693SSebastian Grimberg } 617dce49693SSebastian Grimberg } else if (rstr_type == CEED_RESTRICTION_CURL_ORIENTED) { 618dce49693SSebastian Grimberg switch (mem_type) { 619dce49693SSebastian Grimberg case CEED_MEM_HOST: { 620f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetHostCeedInt8Array(curl_orients, copy_mode, 3 * size, &impl->h_curl_orients_owned, &impl->h_curl_orients_borrowed, 621f5d1e504SJeremy L Thompson &impl->h_curl_orients)); 622a267acd1SJeremy L Thompson CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_curl_orients_owned, 3 * size * sizeof(CeedInt8))); 623f5d1e504SJeremy L Thompson CeedCallCuda(ceed, 624f5d1e504SJeremy L Thompson cudaMemcpy((CeedInt8 *)impl->d_curl_orients_owned, impl->h_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyHostToDevice)); 625a267acd1SJeremy L Thompson impl->d_curl_orients = impl->d_curl_orients_owned; 626dce49693SSebastian Grimberg } break; 627dce49693SSebastian Grimberg case CEED_MEM_DEVICE: { 628f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetDeviceCeedInt8Array_Cuda(ceed, curl_orients, copy_mode, 3 * size, &impl->d_curl_orients_owned, 629f5d1e504SJeremy L Thompson &impl->d_curl_orients_borrowed, (const CeedInt8 **)&impl->d_curl_orients)); 630a267acd1SJeremy L Thompson CeedCallBackend(CeedMalloc(3 * size, &impl->h_curl_orients_owned)); 631f5d1e504SJeremy L Thompson CeedCallCuda(ceed, 632f5d1e504SJeremy L Thompson cudaMemcpy((CeedInt8 *)impl->h_curl_orients_owned, impl->d_curl_orients, 3 * size * sizeof(CeedInt8), cudaMemcpyDeviceToHost)); 633a267acd1SJeremy L Thompson impl->h_curl_orients = impl->h_curl_orients_owned; 634dce49693SSebastian Grimberg } break; 635dce49693SSebastian Grimberg } 636dce49693SSebastian Grimberg } 637dce49693SSebastian Grimberg } 638ca735530SJeremy L Thompson 639ff1e7120SSebastian Grimberg // Register backend functions 640dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Apply", CeedElemRestrictionApply_Cuda)); 641dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnsigned", CeedElemRestrictionApplyUnsigned_Cuda)); 642dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "ApplyUnoriented", CeedElemRestrictionApplyUnoriented_Cuda)); 643dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOffsets", CeedElemRestrictionGetOffsets_Cuda)); 644dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetOrientations", CeedElemRestrictionGetOrientations_Cuda)); 645dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetCurlOrientations", CeedElemRestrictionGetCurlOrientations_Cuda)); 646*b20a4af9SJeremy L Thompson if (rstr_type == CEED_RESTRICTION_POINTS) { 647*b20a4af9SJeremy L Thompson CeedCallBackend( 648*b20a4af9SJeremy L Thompson CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "GetAtPointsElementOffset", CeedElemRestrictionGetAtPointsElementOffset_Cuda)); 649*b20a4af9SJeremy L Thompson } 650dce49693SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "ElemRestriction", rstr, "Destroy", CeedElemRestrictionDestroy_Cuda)); 651ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 652ff1e7120SSebastian Grimberg } 653ff1e7120SSebastian Grimberg 654ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 655