1d275d636SJeremy L Thompson // Copyright (c) 2017-2025, 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 <cuda_runtime.h> 11ff1e7120SSebastian Grimberg #include <math.h> 12ff1e7120SSebastian Grimberg #include <stdbool.h> 13ff1e7120SSebastian Grimberg #include <string.h> 14ff1e7120SSebastian Grimberg 15ff1e7120SSebastian Grimberg #include "../cuda/ceed-cuda-common.h" 16ff1e7120SSebastian Grimberg #include "ceed-cuda-ref.h" 17ff1e7120SSebastian Grimberg 18ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 19ff1e7120SSebastian Grimberg // Check if host/device sync is needed 20ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 21ff1e7120SSebastian Grimberg static inline int CeedVectorNeedSync_Cuda(const CeedVector vec, CeedMemType mem_type, bool *need_sync) { 22ff1e7120SSebastian Grimberg bool has_valid_array = false; 23ca735530SJeremy L Thompson CeedVector_Cuda *impl; 24ca735530SJeremy L Thompson 25ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 26ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorHasValidArray(vec, &has_valid_array)); 27ff1e7120SSebastian Grimberg switch (mem_type) { 28ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 29ff1e7120SSebastian Grimberg *need_sync = has_valid_array && !impl->h_array; 30ff1e7120SSebastian Grimberg break; 31ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 32ff1e7120SSebastian Grimberg *need_sync = has_valid_array && !impl->d_array; 33ff1e7120SSebastian Grimberg break; 34ff1e7120SSebastian Grimberg } 35ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 36ff1e7120SSebastian Grimberg } 37ff1e7120SSebastian Grimberg 38ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 39ff1e7120SSebastian Grimberg // Sync host to device 40ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 41ff1e7120SSebastian Grimberg static inline int CeedVectorSyncH2D_Cuda(const CeedVector vec) { 42ca735530SJeremy L Thompson CeedSize length; 43672b0f2aSSebastian Grimberg size_t bytes; 44ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 45ca735530SJeremy L Thompson 46ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 47ff1e7120SSebastian Grimberg 486e536b99SJeremy L Thompson CeedCheck(impl->h_array, CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "No valid host data to sync to device"); 49ff1e7120SSebastian Grimberg 50ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 51672b0f2aSSebastian Grimberg bytes = length * sizeof(CeedScalar); 52ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) { 53ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_borrowed; 54ff1e7120SSebastian Grimberg } else if (impl->d_array_owned) { 55ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 56ff1e7120SSebastian Grimberg } else { 579bc66399SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaMalloc((void **)&impl->d_array_owned, bytes)); 58ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 59ff1e7120SSebastian Grimberg } 609bc66399SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaMemcpy(impl->d_array, impl->h_array, bytes, cudaMemcpyHostToDevice)); 61ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 62ff1e7120SSebastian Grimberg } 63ff1e7120SSebastian Grimberg 64ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 65ff1e7120SSebastian Grimberg // Sync device to host 66ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 67ff1e7120SSebastian Grimberg static inline int CeedVectorSyncD2H_Cuda(const CeedVector vec) { 68ca735530SJeremy L Thompson CeedSize length; 69ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 70ca735530SJeremy L Thompson 71ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 72ff1e7120SSebastian Grimberg 739bc66399SJeremy L Thompson CeedCheck(impl->d_array, CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "No valid device data to sync to host"); 74ff1e7120SSebastian Grimberg 75ff1e7120SSebastian Grimberg if (impl->h_array_borrowed) { 76ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_borrowed; 77ff1e7120SSebastian Grimberg } else if (impl->h_array_owned) { 78ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 79ff1e7120SSebastian Grimberg } else { 80ff1e7120SSebastian Grimberg CeedSize length; 81ca735530SJeremy L Thompson 82ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 83ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(length, &impl->h_array_owned)); 84ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 85ff1e7120SSebastian Grimberg } 86ff1e7120SSebastian Grimberg 87ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 88ff1e7120SSebastian Grimberg size_t bytes = length * sizeof(CeedScalar); 89ff1e7120SSebastian Grimberg 909bc66399SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaMemcpy(impl->h_array, impl->d_array, bytes, cudaMemcpyDeviceToHost)); 91ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 92ff1e7120SSebastian Grimberg } 93ff1e7120SSebastian Grimberg 94ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 95ff1e7120SSebastian Grimberg // Sync arrays 96ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 97ff1e7120SSebastian Grimberg static int CeedVectorSyncArray_Cuda(const CeedVector vec, CeedMemType mem_type) { 98ff1e7120SSebastian Grimberg bool need_sync = false; 99ca735530SJeremy L Thompson 100ca735530SJeremy L Thompson // Check whether device/host sync is needed 101ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorNeedSync_Cuda(vec, mem_type, &need_sync)); 102ff1e7120SSebastian Grimberg if (!need_sync) return CEED_ERROR_SUCCESS; 103ff1e7120SSebastian Grimberg 104ff1e7120SSebastian Grimberg switch (mem_type) { 105ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 106ff1e7120SSebastian Grimberg return CeedVectorSyncD2H_Cuda(vec); 107ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 108ff1e7120SSebastian Grimberg return CeedVectorSyncH2D_Cuda(vec); 109ff1e7120SSebastian Grimberg } 110ff1e7120SSebastian Grimberg return CEED_ERROR_UNSUPPORTED; 111ff1e7120SSebastian Grimberg } 112ff1e7120SSebastian Grimberg 113ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 114ff1e7120SSebastian Grimberg // Set all pointers as invalid 115ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 116ff1e7120SSebastian Grimberg static inline int CeedVectorSetAllInvalid_Cuda(const CeedVector vec) { 117ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 118ff1e7120SSebastian Grimberg 119ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 120ff1e7120SSebastian Grimberg impl->h_array = NULL; 121ff1e7120SSebastian Grimberg impl->d_array = NULL; 122ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 123ff1e7120SSebastian Grimberg } 124ff1e7120SSebastian Grimberg 125ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 126ff1e7120SSebastian Grimberg // Check if CeedVector has any valid pointer 127ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 128ff1e7120SSebastian Grimberg static inline int CeedVectorHasValidArray_Cuda(const CeedVector vec, bool *has_valid_array) { 129ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 130ca735530SJeremy L Thompson 131ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 1321c66c397SJeremy L Thompson *has_valid_array = impl->h_array || impl->d_array; 133ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 134ff1e7120SSebastian Grimberg } 135ff1e7120SSebastian Grimberg 136ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 137ff1e7120SSebastian Grimberg // Check if has array of given type 138ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 139ff1e7120SSebastian Grimberg static inline int CeedVectorHasArrayOfType_Cuda(const CeedVector vec, CeedMemType mem_type, bool *has_array_of_type) { 140ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 141ff1e7120SSebastian Grimberg 142ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 143ff1e7120SSebastian Grimberg switch (mem_type) { 144ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 1451c66c397SJeremy L Thompson *has_array_of_type = impl->h_array_borrowed || impl->h_array_owned; 146ff1e7120SSebastian Grimberg break; 147ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 1481c66c397SJeremy L Thompson *has_array_of_type = impl->d_array_borrowed || impl->d_array_owned; 149ff1e7120SSebastian Grimberg break; 150ff1e7120SSebastian Grimberg } 151ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 152ff1e7120SSebastian Grimberg } 153ff1e7120SSebastian Grimberg 154ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 155ff1e7120SSebastian Grimberg // Check if has borrowed array of given type 156ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 157ff1e7120SSebastian Grimberg static inline int CeedVectorHasBorrowedArrayOfType_Cuda(const CeedVector vec, CeedMemType mem_type, bool *has_borrowed_array_of_type) { 158ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 159ff1e7120SSebastian Grimberg 160ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 161ff1e7120SSebastian Grimberg switch (mem_type) { 162ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 1631c66c397SJeremy L Thompson *has_borrowed_array_of_type = impl->h_array_borrowed; 164ff1e7120SSebastian Grimberg break; 165ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 1661c66c397SJeremy L Thompson *has_borrowed_array_of_type = impl->d_array_borrowed; 167ff1e7120SSebastian Grimberg break; 168ff1e7120SSebastian Grimberg } 169ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 170ff1e7120SSebastian Grimberg } 171ff1e7120SSebastian Grimberg 172ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 173ff1e7120SSebastian Grimberg // Set array from host 174ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 175ff1e7120SSebastian Grimberg static int CeedVectorSetArrayHost_Cuda(const CeedVector vec, const CeedCopyMode copy_mode, CeedScalar *array) { 176a267acd1SJeremy L Thompson CeedSize length; 177ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 178ff1e7120SSebastian Grimberg 179ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 180ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 181a267acd1SJeremy L Thompson 182f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetHostCeedScalarArray(array, copy_mode, length, (const CeedScalar **)&impl->h_array_owned, 183f5d1e504SJeremy L Thompson (const CeedScalar **)&impl->h_array_borrowed, (const CeedScalar **)&impl->h_array)); 184ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 185ff1e7120SSebastian Grimberg } 186ff1e7120SSebastian Grimberg 187ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 188ff1e7120SSebastian Grimberg // Set array from device 189ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 190ff1e7120SSebastian Grimberg static int CeedVectorSetArrayDevice_Cuda(const CeedVector vec, const CeedCopyMode copy_mode, CeedScalar *array) { 191a267acd1SJeremy L Thompson CeedSize length; 192ff1e7120SSebastian Grimberg Ceed ceed; 193ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 194ff1e7120SSebastian Grimberg 195ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 196ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 197ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 198a267acd1SJeremy L Thompson 199f5d1e504SJeremy L Thompson CeedCallBackend(CeedSetDeviceCeedScalarArray_Cuda(ceed, array, copy_mode, length, (const CeedScalar **)&impl->d_array_owned, 200f5d1e504SJeremy L Thompson (const CeedScalar **)&impl->d_array_borrowed, (const CeedScalar **)&impl->d_array)); 2019bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 202ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 203ff1e7120SSebastian Grimberg } 204ff1e7120SSebastian Grimberg 205ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 206ff1e7120SSebastian Grimberg // Set the array used by a vector, 207ff1e7120SSebastian Grimberg // freeing any previously allocated array if applicable 208ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 209ff1e7120SSebastian Grimberg static int CeedVectorSetArray_Cuda(const CeedVector vec, const CeedMemType mem_type, const CeedCopyMode copy_mode, CeedScalar *array) { 210ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 211ff1e7120SSebastian Grimberg 212ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 213ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetAllInvalid_Cuda(vec)); 214ff1e7120SSebastian Grimberg switch (mem_type) { 215ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 216ff1e7120SSebastian Grimberg return CeedVectorSetArrayHost_Cuda(vec, copy_mode, array); 217ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 218ff1e7120SSebastian Grimberg return CeedVectorSetArrayDevice_Cuda(vec, copy_mode, array); 219ff1e7120SSebastian Grimberg } 220ff1e7120SSebastian Grimberg return CEED_ERROR_UNSUPPORTED; 221ff1e7120SSebastian Grimberg } 222ff1e7120SSebastian Grimberg 223ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 224f1c2287bSJeremy L Thompson // Copy host array to value strided 225f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 226f1c2287bSJeremy L Thompson static int CeedHostCopyStrided_Cuda(CeedScalar *h_array, CeedSize start, CeedSize step, CeedSize length, CeedScalar *h_copy_array) { 227f1c2287bSJeremy L Thompson for (CeedSize i = start; i < length; i += step) h_copy_array[i] = h_array[i]; 228f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 229f1c2287bSJeremy L Thompson } 230f1c2287bSJeremy L Thompson 231f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 232f1c2287bSJeremy L Thompson // Copy device array to value strided (impl in .cu file) 233f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 234f1c2287bSJeremy L Thompson int CeedDeviceCopyStrided_Cuda(CeedScalar *d_array, CeedSize start, CeedSize step, CeedSize length, CeedScalar *d_copy_array); 235f1c2287bSJeremy L Thompson 236f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 237f1c2287bSJeremy L Thompson // Copy a vector to a value strided 238f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 239f1c2287bSJeremy L Thompson static int CeedVectorCopyStrided_Cuda(CeedVector vec, CeedSize start, CeedSize step, CeedVector vec_copy) { 240f1c2287bSJeremy L Thompson CeedSize length; 241f1c2287bSJeremy L Thompson CeedVector_Cuda *impl; 242f1c2287bSJeremy L Thompson 243f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 244a7efc114SJeremy L Thompson { 245a7efc114SJeremy L Thompson CeedSize length_vec, length_copy; 246a7efc114SJeremy L Thompson 2475a5594ffSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length_vec)); 2485a5594ffSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec_copy, &length_copy)); 249a7efc114SJeremy L Thompson length = length_vec < length_copy ? length_vec : length_copy; 250a7efc114SJeremy L Thompson } 251f1c2287bSJeremy L Thompson // Set value for synced device/host array 252f1c2287bSJeremy L Thompson if (impl->d_array) { 253f1c2287bSJeremy L Thompson CeedScalar *copy_array; 254f1c2287bSJeremy L Thompson 255f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetArray(vec_copy, CEED_MEM_DEVICE, ©_array)); 256e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 257e84c3ebcSJeremy L Thompson cublasHandle_t handle; 258e84c3ebcSJeremy L Thompson Ceed ceed; 259e84c3ebcSJeremy L Thompson 260e84c3ebcSJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 261e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(ceed, &handle)); 262e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 263e84c3ebcSJeremy L Thompson CeedCallCublas(ceed, cublasScopy_64(handle, (int64_t)length, impl->d_array + start, (int64_t)step, copy_array + start, (int64_t)step)); 264e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 265e84c3ebcSJeremy L Thompson CeedCallCublas(ceed, cublasDcopy_64(handle, (int64_t)length, impl->d_array + start, (int64_t)step, copy_array + start, (int64_t)step)); 266e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 267e84c3ebcSJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 268e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 269f1c2287bSJeremy L Thompson CeedCallBackend(CeedDeviceCopyStrided_Cuda(impl->d_array, start, step, length, copy_array)); 270e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 271f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(vec_copy, ©_array)); 272e84c3ebcSJeremy L Thompson impl->h_array = NULL; 273f1c2287bSJeremy L Thompson } else if (impl->h_array) { 274f1c2287bSJeremy L Thompson CeedScalar *copy_array; 275f1c2287bSJeremy L Thompson 276f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetArray(vec_copy, CEED_MEM_HOST, ©_array)); 277f1c2287bSJeremy L Thompson CeedCallBackend(CeedHostCopyStrided_Cuda(impl->h_array, start, step, length, copy_array)); 278f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(vec_copy, ©_array)); 279e84c3ebcSJeremy L Thompson impl->d_array = NULL; 280f1c2287bSJeremy L Thompson } else { 281f1c2287bSJeremy L Thompson return CeedError(CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "CeedVector must have valid data set"); 282f1c2287bSJeremy L Thompson } 283f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 284f1c2287bSJeremy L Thompson } 285f1c2287bSJeremy L Thompson 286f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 287ff1e7120SSebastian Grimberg // Set host array to value 288ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 289f7c1b517Snbeams static int CeedHostSetValue_Cuda(CeedScalar *h_array, CeedSize length, CeedScalar val) { 290f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) h_array[i] = val; 291ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 292ff1e7120SSebastian Grimberg } 293ff1e7120SSebastian Grimberg 294ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 295ff1e7120SSebastian Grimberg // Set device array to value (impl in .cu file) 296ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 297f7c1b517Snbeams int CeedDeviceSetValue_Cuda(CeedScalar *d_array, CeedSize length, CeedScalar val); 298ff1e7120SSebastian Grimberg 299ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 300f7c1b517Snbeams // Set a vector to a value 301ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 302ff1e7120SSebastian Grimberg static int CeedVectorSetValue_Cuda(CeedVector vec, CeedScalar val) { 303ff1e7120SSebastian Grimberg CeedSize length; 304ca735530SJeremy L Thompson CeedVector_Cuda *impl; 305ff1e7120SSebastian Grimberg 306ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 307ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 308ff1e7120SSebastian Grimberg // Set value for synced device/host array 309ff1e7120SSebastian Grimberg if (!impl->d_array && !impl->h_array) { 310ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) { 311ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_borrowed; 312ff1e7120SSebastian Grimberg } else if (impl->h_array_borrowed) { 313ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_borrowed; 314ff1e7120SSebastian Grimberg } else if (impl->d_array_owned) { 315ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 316ff1e7120SSebastian Grimberg } else if (impl->h_array_owned) { 317ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 318ff1e7120SSebastian Grimberg } else { 319ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, CEED_MEM_DEVICE, CEED_COPY_VALUES, NULL)); 320ff1e7120SSebastian Grimberg } 321ff1e7120SSebastian Grimberg } 322ff1e7120SSebastian Grimberg if (impl->d_array) { 323124cc107SJeremy L Thompson if (val == 0) { 324124cc107SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaMemset(impl->d_array, 0, length * sizeof(CeedScalar))); 325124cc107SJeremy L Thompson } else { 326ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceSetValue_Cuda(impl->d_array, length, val)); 327ff1e7120SSebastian Grimberg } 328124cc107SJeremy L Thompson impl->h_array = NULL; 329124cc107SJeremy L Thompson } else if (impl->h_array) { 330ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostSetValue_Cuda(impl->h_array, length, val)); 331ff1e7120SSebastian Grimberg impl->d_array = NULL; 332ff1e7120SSebastian Grimberg } 333ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 334ff1e7120SSebastian Grimberg } 335ff1e7120SSebastian Grimberg 336ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 337f1c2287bSJeremy L Thompson // Set host array to value strided 338f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 339*ff90b007SJeremy L Thompson static int CeedHostSetValueStrided_Cuda(CeedScalar *h_array, CeedSize start, CeedSize stop, CeedSize step, CeedSize length, CeedScalar val) { 340*ff90b007SJeremy L Thompson for (CeedSize i = start; i <= stop; i += step) h_array[i] = val; 341f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 342f1c2287bSJeremy L Thompson } 343f1c2287bSJeremy L Thompson 344f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 345f1c2287bSJeremy L Thompson // Set device array to value strided (impl in .cu file) 346f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 347*ff90b007SJeremy L Thompson int CeedDeviceSetValueStrided_Cuda(CeedScalar *d_array, CeedSize start, CeedSize stop, CeedSize step, CeedSize length, CeedScalar val); 348f1c2287bSJeremy L Thompson 349f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 350f1c2287bSJeremy L Thompson // Set a vector to a value strided 351f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 352*ff90b007SJeremy L Thompson static int CeedVectorSetValueStrided_Cuda(CeedVector vec, CeedSize start, CeedSize stop, CeedSize step, CeedScalar val) { 353f1c2287bSJeremy L Thompson CeedSize length; 354f1c2287bSJeremy L Thompson CeedVector_Cuda *impl; 355f1c2287bSJeremy L Thompson 356f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 357f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 358f1c2287bSJeremy L Thompson // Set value for synced device/host array 359*ff90b007SJeremy L Thompson if (stop == -1) stop = length; 360f1c2287bSJeremy L Thompson if (impl->d_array) { 361*ff90b007SJeremy L Thompson CeedCallBackend(CeedDeviceSetValueStrided_Cuda(impl->d_array, start, stop, step, length, val)); 362f1c2287bSJeremy L Thompson impl->h_array = NULL; 363f1c2287bSJeremy L Thompson } else if (impl->h_array) { 364*ff90b007SJeremy L Thompson CeedCallBackend(CeedHostSetValueStrided_Cuda(impl->h_array, start, stop, step, length, val)); 365f1c2287bSJeremy L Thompson impl->d_array = NULL; 366f1c2287bSJeremy L Thompson } else { 367f1c2287bSJeremy L Thompson return CeedError(CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "CeedVector must have valid data set"); 368f1c2287bSJeremy L Thompson } 369f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 370f1c2287bSJeremy L Thompson } 371f1c2287bSJeremy L Thompson 372f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 373ff1e7120SSebastian Grimberg // Vector Take Array 374ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 375ff1e7120SSebastian Grimberg static int CeedVectorTakeArray_Cuda(CeedVector vec, CeedMemType mem_type, CeedScalar **array) { 376ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 377ff1e7120SSebastian Grimberg 378ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 379ff1e7120SSebastian Grimberg // Sync array to requested mem_type 380ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 381ff1e7120SSebastian Grimberg // Update pointer 382ff1e7120SSebastian Grimberg switch (mem_type) { 383ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 384ff1e7120SSebastian Grimberg (*array) = impl->h_array_borrowed; 385ff1e7120SSebastian Grimberg impl->h_array_borrowed = NULL; 386ff1e7120SSebastian Grimberg impl->h_array = NULL; 387ff1e7120SSebastian Grimberg break; 388ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 389ff1e7120SSebastian Grimberg (*array) = impl->d_array_borrowed; 390ff1e7120SSebastian Grimberg impl->d_array_borrowed = NULL; 391ff1e7120SSebastian Grimberg impl->d_array = NULL; 392ff1e7120SSebastian Grimberg break; 393ff1e7120SSebastian Grimberg } 394ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 395ff1e7120SSebastian Grimberg } 396ff1e7120SSebastian Grimberg 397ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 398ff1e7120SSebastian Grimberg // Core logic for array syncronization for GetArray. 399ff1e7120SSebastian Grimberg // If a different memory type is most up to date, this will perform a copy 400ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 401ff1e7120SSebastian Grimberg static int CeedVectorGetArrayCore_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 402ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 403ff1e7120SSebastian Grimberg 404ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 405ff1e7120SSebastian Grimberg // Sync array to requested mem_type 406ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 407ff1e7120SSebastian Grimberg // Update pointer 408ff1e7120SSebastian Grimberg switch (mem_type) { 409ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 410ff1e7120SSebastian Grimberg *array = impl->h_array; 411ff1e7120SSebastian Grimberg break; 412ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 413ff1e7120SSebastian Grimberg *array = impl->d_array; 414ff1e7120SSebastian Grimberg break; 415ff1e7120SSebastian Grimberg } 416ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 417ff1e7120SSebastian Grimberg } 418ff1e7120SSebastian Grimberg 419ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 420ff1e7120SSebastian Grimberg // Get read-only access to a vector via the specified mem_type 421ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 422ff1e7120SSebastian Grimberg static int CeedVectorGetArrayRead_Cuda(const CeedVector vec, const CeedMemType mem_type, const CeedScalar **array) { 423ff1e7120SSebastian Grimberg return CeedVectorGetArrayCore_Cuda(vec, mem_type, (CeedScalar **)array); 424ff1e7120SSebastian Grimberg } 425ff1e7120SSebastian Grimberg 426ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 427ff1e7120SSebastian Grimberg // Get read/write access to a vector via the specified mem_type 428ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 429ff1e7120SSebastian Grimberg static int CeedVectorGetArray_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 430ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 431ca735530SJeremy L Thompson 432ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 433ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayCore_Cuda(vec, mem_type, array)); 434ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetAllInvalid_Cuda(vec)); 435ff1e7120SSebastian Grimberg switch (mem_type) { 436ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 437ff1e7120SSebastian Grimberg impl->h_array = *array; 438ff1e7120SSebastian Grimberg break; 439ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 440ff1e7120SSebastian Grimberg impl->d_array = *array; 441ff1e7120SSebastian Grimberg break; 442ff1e7120SSebastian Grimberg } 443ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 444ff1e7120SSebastian Grimberg } 445ff1e7120SSebastian Grimberg 446ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 447ff1e7120SSebastian Grimberg // Get write access to a vector via the specified mem_type 448ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 449ff1e7120SSebastian Grimberg static int CeedVectorGetArrayWrite_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 450ff1e7120SSebastian Grimberg bool has_array_of_type = true; 451ca735530SJeremy L Thompson CeedVector_Cuda *impl; 452ca735530SJeremy L Thompson 453ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 454ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorHasArrayOfType_Cuda(vec, mem_type, &has_array_of_type)); 455ff1e7120SSebastian Grimberg if (!has_array_of_type) { 456ff1e7120SSebastian Grimberg // Allocate if array is not yet allocated 457ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, mem_type, CEED_COPY_VALUES, NULL)); 458ff1e7120SSebastian Grimberg } else { 459ff1e7120SSebastian Grimberg // Select dirty array 460ff1e7120SSebastian Grimberg switch (mem_type) { 461ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 462ff1e7120SSebastian Grimberg if (impl->h_array_borrowed) impl->h_array = impl->h_array_borrowed; 463ff1e7120SSebastian Grimberg else impl->h_array = impl->h_array_owned; 464ff1e7120SSebastian Grimberg break; 465ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 466ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) impl->d_array = impl->d_array_borrowed; 467ff1e7120SSebastian Grimberg else impl->d_array = impl->d_array_owned; 468ff1e7120SSebastian Grimberg } 469ff1e7120SSebastian Grimberg } 470ff1e7120SSebastian Grimberg return CeedVectorGetArray_Cuda(vec, mem_type, array); 471ff1e7120SSebastian Grimberg } 472ff1e7120SSebastian Grimberg 473ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 474ff1e7120SSebastian Grimberg // Get the norm of a CeedVector 475ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 476ff1e7120SSebastian Grimberg static int CeedVectorNorm_Cuda(CeedVector vec, CeedNormType type, CeedScalar *norm) { 477ff1e7120SSebastian Grimberg Ceed ceed; 478ca735530SJeremy L Thompson CeedSize length; 479e84c3ebcSJeremy L Thompson #if (CUDA_VERSION < 12000) 480672b0f2aSSebastian Grimberg CeedSize num_calls; 481e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 482ca735530SJeremy L Thompson const CeedScalar *d_array; 483ca735530SJeremy L Thompson CeedVector_Cuda *impl; 484672b0f2aSSebastian Grimberg cublasHandle_t handle; 485ca735530SJeremy L Thompson 486ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 487ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 488ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 489ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetCublasHandle_Cuda(ceed, &handle)); 490ff1e7120SSebastian Grimberg 491e84c3ebcSJeremy L Thompson #if (CUDA_VERSION < 12000) 492f7c1b517Snbeams // With CUDA 12, we can use the 64-bit integer interface. Prior to that, 493f7c1b517Snbeams // we need to check if the vector is too long to handle with int32, 494b2165e7aSSebastian Grimberg // and if so, divide it into subsections for repeated cuBLAS calls. 495672b0f2aSSebastian Grimberg num_calls = length / INT_MAX; 496f7c1b517Snbeams if (length % INT_MAX > 0) num_calls += 1; 497e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 498f7c1b517Snbeams 499ff1e7120SSebastian Grimberg // Compute norm 500ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &d_array)); 501ff1e7120SSebastian Grimberg switch (type) { 502ff1e7120SSebastian Grimberg case CEED_NORM_1: { 503f6f49adbSnbeams *norm = 0.0; 504e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 505e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) // We have CUDA 12, and can use 64-bit integers 506f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 507e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 508f7c1b517Snbeams float sub_norm = 0.0; 509f7c1b517Snbeams float *d_array_start; 510ca735530SJeremy L Thompson 511f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 512f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 513f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 514f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 515ca735530SJeremy L Thompson 516f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 517f7c1b517Snbeams *norm += sub_norm; 518f7c1b517Snbeams } 519e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 520e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 521e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 522f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 523e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 524f7c1b517Snbeams double sub_norm = 0.0; 525f7c1b517Snbeams double *d_array_start; 526ca735530SJeremy L Thompson 527f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 528f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 529f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 530f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 531ca735530SJeremy L Thompson 532f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 533f7c1b517Snbeams *norm += sub_norm; 534f7c1b517Snbeams } 535e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 536e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 537ff1e7120SSebastian Grimberg break; 538ff1e7120SSebastian Grimberg } 539ff1e7120SSebastian Grimberg case CEED_NORM_2: { 540e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 541e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 542f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 543e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 544f7c1b517Snbeams float sub_norm = 0.0, norm_sum = 0.0; 545f7c1b517Snbeams float *d_array_start; 546ca735530SJeremy L Thompson 547f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 548f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 549f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 550f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 551ca735530SJeremy L Thompson 552f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 553f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 554f7c1b517Snbeams } 555f7c1b517Snbeams *norm = sqrt(norm_sum); 556e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 557e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 558e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 559f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 560e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 561f7c1b517Snbeams double sub_norm = 0.0, norm_sum = 0.0; 562f7c1b517Snbeams double *d_array_start; 563ca735530SJeremy L Thompson 564f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 565f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 566f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 567f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 568ca735530SJeremy L Thompson 569f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 570f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 571f7c1b517Snbeams } 572f7c1b517Snbeams *norm = sqrt(norm_sum); 573e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 574e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 575ff1e7120SSebastian Grimberg break; 576ff1e7120SSebastian Grimberg } 577ff1e7120SSebastian Grimberg case CEED_NORM_MAX: { 578e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 579e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 580ca735530SJeremy L Thompson int64_t index; 581ca735530SJeremy L Thompson CeedScalar norm_no_abs; 582ca735530SJeremy L Thompson 583ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax_64(handle, (int64_t)length, (float *)d_array, 1, &index)); 584ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 585ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 586e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 587ca735530SJeremy L Thompson CeedInt index; 588f7c1b517Snbeams float sub_max = 0.0, current_max = 0.0; 589f7c1b517Snbeams float *d_array_start; 590ca735530SJeremy L Thompson 591f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 592f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 593f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 594f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 595ca735530SJeremy L Thompson 596ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &index)); 597ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 598f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 599f7c1b517Snbeams } 600f7c1b517Snbeams *norm = current_max; 601e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 602e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 603e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 604ca735530SJeremy L Thompson int64_t index; 605ca735530SJeremy L Thompson CeedScalar norm_no_abs; 606ca735530SJeremy L Thompson 607ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax_64(handle, (int64_t)length, (double *)d_array, 1, &index)); 608ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 609ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 610e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 611ca735530SJeremy L Thompson CeedInt index; 612f7c1b517Snbeams double sub_max = 0.0, current_max = 0.0; 613f7c1b517Snbeams double *d_array_start; 614ca735530SJeremy L Thompson 615f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 616f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 617f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 618f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 619ca735530SJeremy L Thompson 620ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &index)); 621ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 622f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 623f7c1b517Snbeams } 624f7c1b517Snbeams *norm = current_max; 625e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 626e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 627ff1e7120SSebastian Grimberg break; 628ff1e7120SSebastian Grimberg } 629ff1e7120SSebastian Grimberg } 630ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(vec, &d_array)); 6319bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 632ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 633ff1e7120SSebastian Grimberg } 634ff1e7120SSebastian Grimberg 635ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 636ff1e7120SSebastian Grimberg // Take reciprocal of a vector on host 637ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 638f7c1b517Snbeams static int CeedHostReciprocal_Cuda(CeedScalar *h_array, CeedSize length) { 639f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) { 640ff1e7120SSebastian Grimberg if (fabs(h_array[i]) > CEED_EPSILON) h_array[i] = 1. / h_array[i]; 641ff1e7120SSebastian Grimberg } 642ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 643ff1e7120SSebastian Grimberg } 644ff1e7120SSebastian Grimberg 645ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 646ff1e7120SSebastian Grimberg // Take reciprocal of a vector on device (impl in .cu file) 647ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 648f7c1b517Snbeams int CeedDeviceReciprocal_Cuda(CeedScalar *d_array, CeedSize length); 649ff1e7120SSebastian Grimberg 650ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 651ff1e7120SSebastian Grimberg // Take reciprocal of a vector 652ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 653ff1e7120SSebastian Grimberg static int CeedVectorReciprocal_Cuda(CeedVector vec) { 654ff1e7120SSebastian Grimberg CeedSize length; 655ca735530SJeremy L Thompson CeedVector_Cuda *impl; 656ff1e7120SSebastian Grimberg 657ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 658ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 659ff1e7120SSebastian Grimberg // Set value for synced device/host array 660ff1e7120SSebastian Grimberg if (impl->d_array) CeedCallBackend(CeedDeviceReciprocal_Cuda(impl->d_array, length)); 661ff1e7120SSebastian Grimberg if (impl->h_array) CeedCallBackend(CeedHostReciprocal_Cuda(impl->h_array, length)); 662ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 663ff1e7120SSebastian Grimberg } 664ff1e7120SSebastian Grimberg 665ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 666ff1e7120SSebastian Grimberg // Compute x = alpha x on the host 667ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 668f7c1b517Snbeams static int CeedHostScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length) { 669f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) x_array[i] *= alpha; 670ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 671ff1e7120SSebastian Grimberg } 672ff1e7120SSebastian Grimberg 673ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 674ff1e7120SSebastian Grimberg // Compute x = alpha x on device (impl in .cu file) 675ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 676f7c1b517Snbeams int CeedDeviceScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length); 677ff1e7120SSebastian Grimberg 678ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 679ff1e7120SSebastian Grimberg // Compute x = alpha x 680ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 681ff1e7120SSebastian Grimberg static int CeedVectorScale_Cuda(CeedVector x, CeedScalar alpha) { 682ff1e7120SSebastian Grimberg CeedSize length; 683e84c3ebcSJeremy L Thompson CeedVector_Cuda *impl; 684ff1e7120SSebastian Grimberg 685e84c3ebcSJeremy L Thompson CeedCallBackend(CeedVectorGetData(x, &impl)); 686ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(x, &length)); 687ff1e7120SSebastian Grimberg // Set value for synced device/host array 688e84c3ebcSJeremy L Thompson if (impl->d_array) { 689e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 690e84c3ebcSJeremy L Thompson cublasHandle_t handle; 691e84c3ebcSJeremy L Thompson 692e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(CeedVectorReturnCeed(x), &handle)); 693e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 694e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(x), cublasSscal_64(handle, (int64_t)length, &alpha, impl->d_array, 1)); 695e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 696e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(x), cublasDscal_64(handle, (int64_t)length, &alpha, impl->d_array, 1)); 697e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 698e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 699e84c3ebcSJeremy L Thompson CeedCallBackend(CeedDeviceScale_Cuda(impl->d_array, alpha, length)); 700e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 701e84c3ebcSJeremy L Thompson impl->h_array = NULL; 702e84c3ebcSJeremy L Thompson } else if (impl->h_array) { 703e84c3ebcSJeremy L Thompson CeedCallBackend(CeedHostScale_Cuda(impl->h_array, alpha, length)); 704e84c3ebcSJeremy L Thompson impl->d_array = NULL; 705e84c3ebcSJeremy L Thompson } 706ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 707ff1e7120SSebastian Grimberg } 708ff1e7120SSebastian Grimberg 709ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 710ff1e7120SSebastian Grimberg // Compute y = alpha x + y on the host 711ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 712f7c1b517Snbeams static int CeedHostAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length) { 713f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) y_array[i] += alpha * x_array[i]; 714ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 715ff1e7120SSebastian Grimberg } 716ff1e7120SSebastian Grimberg 717ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 718ff1e7120SSebastian Grimberg // Compute y = alpha x + y on device (impl in .cu file) 719ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 720f7c1b517Snbeams int CeedDeviceAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length); 721ff1e7120SSebastian Grimberg 722ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 723ff1e7120SSebastian Grimberg // Compute y = alpha x + y 724ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 725ff1e7120SSebastian Grimberg static int CeedVectorAXPY_Cuda(CeedVector y, CeedScalar alpha, CeedVector x) { 726ca735530SJeremy L Thompson CeedSize length; 727ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 728ca735530SJeremy L Thompson 729ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 730ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 731ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 732ff1e7120SSebastian Grimberg // Set value for synced device/host array 733ff1e7120SSebastian Grimberg if (y_impl->d_array) { 734ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 735e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 736e84c3ebcSJeremy L Thompson cublasHandle_t handle; 737e84c3ebcSJeremy L Thompson 738e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(CeedVectorReturnCeed(y), &handle)); 739e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 740e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(y), cublasSaxpy_64(handle, (int64_t)length, &alpha, x_impl->d_array, 1, y_impl->d_array, 1)); 741e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 742e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(y), cublasDaxpy_64(handle, (int64_t)length, &alpha, x_impl->d_array, 1, y_impl->d_array, 1)); 743e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 744e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 745ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPY_Cuda(y_impl->d_array, alpha, x_impl->d_array, length)); 746e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 747e84c3ebcSJeremy L Thompson y_impl->h_array = NULL; 748e84c3ebcSJeremy L Thompson } else if (y_impl->h_array) { 749ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 750ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPY_Cuda(y_impl->h_array, alpha, x_impl->h_array, length)); 751e84c3ebcSJeremy L Thompson y_impl->d_array = NULL; 752ff1e7120SSebastian Grimberg } 753ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 754ff1e7120SSebastian Grimberg } 755ff1e7120SSebastian Grimberg 756ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 757ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on the host 758ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 759f7c1b517Snbeams static int CeedHostAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length) { 760aa67b842SZach Atkins for (CeedSize i = 0; i < length; i++) y_array[i] = alpha * x_array[i] + beta * y_array[i]; 761ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 762ff1e7120SSebastian Grimberg } 763ff1e7120SSebastian Grimberg 764ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 765ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on device (impl in .cu file) 766ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 767f7c1b517Snbeams int CeedDeviceAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length); 768ff1e7120SSebastian Grimberg 769ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 770ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y 771ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 772ff1e7120SSebastian Grimberg static int CeedVectorAXPBY_Cuda(CeedVector y, CeedScalar alpha, CeedScalar beta, CeedVector x) { 773ca735530SJeremy L Thompson CeedSize length; 774ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 775ca735530SJeremy L Thompson 776ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 777ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 778ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 779ff1e7120SSebastian Grimberg // Set value for synced device/host array 780ff1e7120SSebastian Grimberg if (y_impl->d_array) { 781ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 782ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPBY_Cuda(y_impl->d_array, alpha, beta, x_impl->d_array, length)); 783ff1e7120SSebastian Grimberg } 784ff1e7120SSebastian Grimberg if (y_impl->h_array) { 785ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 786ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPBY_Cuda(y_impl->h_array, alpha, beta, x_impl->h_array, length)); 787ff1e7120SSebastian Grimberg } 788ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 789ff1e7120SSebastian Grimberg } 790ff1e7120SSebastian Grimberg 791ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 792ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on the host 793ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 794f7c1b517Snbeams static int CeedHostPointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length) { 795f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) w_array[i] = x_array[i] * y_array[i]; 796ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 797ff1e7120SSebastian Grimberg } 798ff1e7120SSebastian Grimberg 799ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 800ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on device (impl in .cu file) 801ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 802f7c1b517Snbeams int CeedDevicePointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length); 803ff1e7120SSebastian Grimberg 804ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 805ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y 806ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 807ff1e7120SSebastian Grimberg static int CeedVectorPointwiseMult_Cuda(CeedVector w, CeedVector x, CeedVector y) { 808ca735530SJeremy L Thompson CeedSize length; 809ff1e7120SSebastian Grimberg CeedVector_Cuda *w_impl, *x_impl, *y_impl; 810ca735530SJeremy L Thompson 811ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(w, &w_impl)); 812ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 813ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 814ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(w, &length)); 815ff1e7120SSebastian Grimberg // Set value for synced device/host array 816ff1e7120SSebastian Grimberg if (!w_impl->d_array && !w_impl->h_array) { 817ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetValue(w, 0.0)); 818ff1e7120SSebastian Grimberg } 819ff1e7120SSebastian Grimberg if (w_impl->d_array) { 820ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 821ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_DEVICE)); 822ff1e7120SSebastian Grimberg CeedCallBackend(CeedDevicePointwiseMult_Cuda(w_impl->d_array, x_impl->d_array, y_impl->d_array, length)); 823ff1e7120SSebastian Grimberg } 824ff1e7120SSebastian Grimberg if (w_impl->h_array) { 825ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 826ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_HOST)); 827ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostPointwiseMult_Cuda(w_impl->h_array, x_impl->h_array, y_impl->h_array, length)); 828ff1e7120SSebastian Grimberg } 829ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 830ff1e7120SSebastian Grimberg } 831ff1e7120SSebastian Grimberg 832ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 833ff1e7120SSebastian Grimberg // Destroy the vector 834ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 835ff1e7120SSebastian Grimberg static int CeedVectorDestroy_Cuda(const CeedVector vec) { 836ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 837ff1e7120SSebastian Grimberg 838ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 8396e536b99SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaFree(impl->d_array_owned)); 840ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_array_owned)); 841ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 842ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 843ff1e7120SSebastian Grimberg } 844ff1e7120SSebastian Grimberg 845ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 846ff1e7120SSebastian Grimberg // Create a vector of the specified length (does not allocate memory) 847ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 848ff1e7120SSebastian Grimberg int CeedVectorCreate_Cuda(CeedSize n, CeedVector vec) { 849ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 850ff1e7120SSebastian Grimberg Ceed ceed; 851ff1e7120SSebastian Grimberg 852ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 853ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasValidArray", CeedVectorHasValidArray_Cuda)); 854ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasBorrowedArrayOfType", CeedVectorHasBorrowedArrayOfType_Cuda)); 855ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetArray", CeedVectorSetArray_Cuda)); 856ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "TakeArray", CeedVectorTakeArray_Cuda)); 8573e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "CopyStrided", CeedVectorCopyStrided_Cuda)); 8583e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetValue", CeedVectorSetValue_Cuda)); 8593e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetValueStrided", CeedVectorSetValueStrided_Cuda)); 860ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SyncArray", CeedVectorSyncArray_Cuda)); 861ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArray", CeedVectorGetArray_Cuda)); 862ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayRead", CeedVectorGetArrayRead_Cuda)); 863ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayWrite", CeedVectorGetArrayWrite_Cuda)); 864ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Norm", CeedVectorNorm_Cuda)); 865ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Reciprocal", CeedVectorReciprocal_Cuda)); 8663e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Scale", CeedVectorScale_Cuda)); 8673e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPY", CeedVectorAXPY_Cuda)); 8683e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPBY", CeedVectorAXPBY_Cuda)); 869ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "PointwiseMult", CeedVectorPointwiseMult_Cuda)); 870ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Destroy", CeedVectorDestroy_Cuda)); 8719bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 872ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 873ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetData(vec, impl)); 874ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 875ff1e7120SSebastian Grimberg } 876ff1e7120SSebastian Grimberg 877ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 878