1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 //------------------------------------------------------------------------------ 226832a6d73SJeremy L Thompson static int CeedHostCopyStrided_Cuda(CeedScalar *h_array, CeedSize start, CeedSize stop, CeedSize step, CeedScalar *h_copy_array) { 227832a6d73SJeremy L Thompson for (CeedSize i = start; i < stop; 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 //------------------------------------------------------------------------------ 234832a6d73SJeremy L Thompson int CeedDeviceCopyStrided_Cuda(CeedScalar *d_array, CeedSize start, CeedSize stop, CeedSize step, CeedScalar *d_copy_array); 235f1c2287bSJeremy L Thompson 236f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 237f1c2287bSJeremy L Thompson // Copy a vector to a value strided 238f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 239832a6d73SJeremy L Thompson static int CeedVectorCopyStrided_Cuda(CeedVector vec, CeedSize start, CeedSize stop, 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 } 251832a6d73SJeremy L Thompson if (stop == -1) stop = length; 252f1c2287bSJeremy L Thompson // Set value for synced device/host array 253f1c2287bSJeremy L Thompson if (impl->d_array) { 254f1c2287bSJeremy L Thompson CeedScalar *copy_array; 255f1c2287bSJeremy L Thompson 256f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetArray(vec_copy, CEED_MEM_DEVICE, ©_array)); 257e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 258e84c3ebcSJeremy L Thompson cublasHandle_t handle; 259e84c3ebcSJeremy L Thompson Ceed ceed; 260e84c3ebcSJeremy L Thompson 261e84c3ebcSJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 262e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(ceed, &handle)); 263e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 264832a6d73SJeremy L Thompson CeedCallCublas(ceed, cublasScopy_64(handle, (int64_t)(stop - start), impl->d_array + start, (int64_t)step, copy_array + start, (int64_t)step)); 265e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 266832a6d73SJeremy L Thompson CeedCallCublas(ceed, cublasDcopy_64(handle, (int64_t)(stop - start), impl->d_array + start, (int64_t)step, copy_array + start, (int64_t)step)); 267e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 268e84c3ebcSJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 269e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 270832a6d73SJeremy L Thompson CeedCallBackend(CeedDeviceCopyStrided_Cuda(impl->d_array, start, stop, step, copy_array)); 271e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 272f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(vec_copy, ©_array)); 273e84c3ebcSJeremy L Thompson impl->h_array = NULL; 274f1c2287bSJeremy L Thompson } else if (impl->h_array) { 275f1c2287bSJeremy L Thompson CeedScalar *copy_array; 276f1c2287bSJeremy L Thompson 277f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetArray(vec_copy, CEED_MEM_HOST, ©_array)); 278832a6d73SJeremy L Thompson CeedCallBackend(CeedHostCopyStrided_Cuda(impl->h_array, start, stop, step, copy_array)); 279f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorRestoreArray(vec_copy, ©_array)); 280e84c3ebcSJeremy L Thompson impl->d_array = NULL; 281f1c2287bSJeremy L Thompson } else { 282f1c2287bSJeremy L Thompson return CeedError(CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "CeedVector must have valid data set"); 283f1c2287bSJeremy L Thompson } 284f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 285f1c2287bSJeremy L Thompson } 286f1c2287bSJeremy L Thompson 287f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 288ff1e7120SSebastian Grimberg // Set host array to value 289ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 290f7c1b517Snbeams static int CeedHostSetValue_Cuda(CeedScalar *h_array, CeedSize length, CeedScalar val) { 291f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) h_array[i] = val; 292ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 293ff1e7120SSebastian Grimberg } 294ff1e7120SSebastian Grimberg 295ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 296ff1e7120SSebastian Grimberg // Set device array to value (impl in .cu file) 297ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 298f7c1b517Snbeams int CeedDeviceSetValue_Cuda(CeedScalar *d_array, CeedSize length, CeedScalar val); 299ff1e7120SSebastian Grimberg 300ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 301f7c1b517Snbeams // Set a vector to a value 302ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 303ff1e7120SSebastian Grimberg static int CeedVectorSetValue_Cuda(CeedVector vec, CeedScalar val) { 304ff1e7120SSebastian Grimberg CeedSize length; 305ca735530SJeremy L Thompson CeedVector_Cuda *impl; 306ff1e7120SSebastian Grimberg 307ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 308ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 309ff1e7120SSebastian Grimberg // Set value for synced device/host array 310ff1e7120SSebastian Grimberg if (!impl->d_array && !impl->h_array) { 311ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) { 312ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_borrowed; 313ff1e7120SSebastian Grimberg } else if (impl->h_array_borrowed) { 314ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_borrowed; 315ff1e7120SSebastian Grimberg } else if (impl->d_array_owned) { 316ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 317ff1e7120SSebastian Grimberg } else if (impl->h_array_owned) { 318ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 319ff1e7120SSebastian Grimberg } else { 320ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, CEED_MEM_DEVICE, CEED_COPY_VALUES, NULL)); 321ff1e7120SSebastian Grimberg } 322ff1e7120SSebastian Grimberg } 323ff1e7120SSebastian Grimberg if (impl->d_array) { 324124cc107SJeremy L Thompson if (val == 0) { 325124cc107SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaMemset(impl->d_array, 0, length * sizeof(CeedScalar))); 326124cc107SJeremy L Thompson } else { 327ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceSetValue_Cuda(impl->d_array, length, val)); 328ff1e7120SSebastian Grimberg } 329124cc107SJeremy L Thompson impl->h_array = NULL; 330124cc107SJeremy L Thompson } else if (impl->h_array) { 331ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostSetValue_Cuda(impl->h_array, length, val)); 332ff1e7120SSebastian Grimberg impl->d_array = NULL; 333ff1e7120SSebastian Grimberg } 334ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 335ff1e7120SSebastian Grimberg } 336ff1e7120SSebastian Grimberg 337ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 338f1c2287bSJeremy L Thompson // Set host array to value strided 339f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 34014c82621SJeremy L Thompson static int CeedHostSetValueStrided_Cuda(CeedScalar *h_array, CeedSize start, CeedSize stop, CeedSize step, CeedScalar val) { 3412d73a370SJeremy L Thompson for (CeedSize i = start; i < stop; i += step) h_array[i] = val; 342f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 343f1c2287bSJeremy L Thompson } 344f1c2287bSJeremy L Thompson 345f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 346f1c2287bSJeremy L Thompson // Set device array to value strided (impl in .cu file) 347f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 34814c82621SJeremy L Thompson int CeedDeviceSetValueStrided_Cuda(CeedScalar *d_array, CeedSize start, CeedSize stop, CeedSize step, CeedScalar val); 349f1c2287bSJeremy L Thompson 350f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 351f1c2287bSJeremy L Thompson // Set a vector to a value strided 352f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 353ff90b007SJeremy L Thompson static int CeedVectorSetValueStrided_Cuda(CeedVector vec, CeedSize start, CeedSize stop, CeedSize step, CeedScalar val) { 354f1c2287bSJeremy L Thompson CeedSize length; 355f1c2287bSJeremy L Thompson CeedVector_Cuda *impl; 356f1c2287bSJeremy L Thompson 357f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 358f1c2287bSJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 359f1c2287bSJeremy L Thompson // Set value for synced device/host array 360ff90b007SJeremy L Thompson if (stop == -1) stop = length; 361f1c2287bSJeremy L Thompson if (impl->d_array) { 36214c82621SJeremy L Thompson CeedCallBackend(CeedDeviceSetValueStrided_Cuda(impl->d_array, start, stop, step, val)); 363f1c2287bSJeremy L Thompson impl->h_array = NULL; 364f1c2287bSJeremy L Thompson } else if (impl->h_array) { 36514c82621SJeremy L Thompson CeedCallBackend(CeedHostSetValueStrided_Cuda(impl->h_array, start, stop, step, val)); 366f1c2287bSJeremy L Thompson impl->d_array = NULL; 367f1c2287bSJeremy L Thompson } else { 368f1c2287bSJeremy L Thompson return CeedError(CeedVectorReturnCeed(vec), CEED_ERROR_BACKEND, "CeedVector must have valid data set"); 369f1c2287bSJeremy L Thompson } 370f1c2287bSJeremy L Thompson return CEED_ERROR_SUCCESS; 371f1c2287bSJeremy L Thompson } 372f1c2287bSJeremy L Thompson 373f1c2287bSJeremy L Thompson //------------------------------------------------------------------------------ 374ff1e7120SSebastian Grimberg // Vector Take Array 375ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 376ff1e7120SSebastian Grimberg static int CeedVectorTakeArray_Cuda(CeedVector vec, CeedMemType mem_type, CeedScalar **array) { 377ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 378ff1e7120SSebastian Grimberg 379ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 380ff1e7120SSebastian Grimberg // Sync array to requested mem_type 381ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 382ff1e7120SSebastian Grimberg // Update pointer 383ff1e7120SSebastian Grimberg switch (mem_type) { 384ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 385ff1e7120SSebastian Grimberg (*array) = impl->h_array_borrowed; 386ff1e7120SSebastian Grimberg impl->h_array_borrowed = NULL; 387ff1e7120SSebastian Grimberg impl->h_array = NULL; 388ff1e7120SSebastian Grimberg break; 389ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 390ff1e7120SSebastian Grimberg (*array) = impl->d_array_borrowed; 391ff1e7120SSebastian Grimberg impl->d_array_borrowed = NULL; 392ff1e7120SSebastian Grimberg impl->d_array = NULL; 393ff1e7120SSebastian Grimberg break; 394ff1e7120SSebastian Grimberg } 395ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 396ff1e7120SSebastian Grimberg } 397ff1e7120SSebastian Grimberg 398ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 399ff1e7120SSebastian Grimberg // Core logic for array syncronization for GetArray. 400ff1e7120SSebastian Grimberg // If a different memory type is most up to date, this will perform a copy 401ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 402ff1e7120SSebastian Grimberg static int CeedVectorGetArrayCore_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 403ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 404ff1e7120SSebastian Grimberg 405ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 406ff1e7120SSebastian Grimberg // Sync array to requested mem_type 407ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 408ff1e7120SSebastian Grimberg // Update pointer 409ff1e7120SSebastian Grimberg switch (mem_type) { 410ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 411ff1e7120SSebastian Grimberg *array = impl->h_array; 412ff1e7120SSebastian Grimberg break; 413ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 414ff1e7120SSebastian Grimberg *array = impl->d_array; 415ff1e7120SSebastian Grimberg break; 416ff1e7120SSebastian Grimberg } 417ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 418ff1e7120SSebastian Grimberg } 419ff1e7120SSebastian Grimberg 420ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 421ff1e7120SSebastian Grimberg // Get read-only access to a vector via the specified mem_type 422ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 423ff1e7120SSebastian Grimberg static int CeedVectorGetArrayRead_Cuda(const CeedVector vec, const CeedMemType mem_type, const CeedScalar **array) { 424ff1e7120SSebastian Grimberg return CeedVectorGetArrayCore_Cuda(vec, mem_type, (CeedScalar **)array); 425ff1e7120SSebastian Grimberg } 426ff1e7120SSebastian Grimberg 427ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 428ff1e7120SSebastian Grimberg // Get read/write access to a vector via the specified mem_type 429ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 430ff1e7120SSebastian Grimberg static int CeedVectorGetArray_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 431ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 432ca735530SJeremy L Thompson 433ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 434ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayCore_Cuda(vec, mem_type, array)); 435ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetAllInvalid_Cuda(vec)); 436ff1e7120SSebastian Grimberg switch (mem_type) { 437ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 438ff1e7120SSebastian Grimberg impl->h_array = *array; 439ff1e7120SSebastian Grimberg break; 440ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 441ff1e7120SSebastian Grimberg impl->d_array = *array; 442ff1e7120SSebastian Grimberg break; 443ff1e7120SSebastian Grimberg } 444ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 445ff1e7120SSebastian Grimberg } 446ff1e7120SSebastian Grimberg 447ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 448ff1e7120SSebastian Grimberg // Get write access to a vector via the specified mem_type 449ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 450ff1e7120SSebastian Grimberg static int CeedVectorGetArrayWrite_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 451ff1e7120SSebastian Grimberg bool has_array_of_type = true; 452ca735530SJeremy L Thompson CeedVector_Cuda *impl; 453ca735530SJeremy L Thompson 454ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 455ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorHasArrayOfType_Cuda(vec, mem_type, &has_array_of_type)); 456ff1e7120SSebastian Grimberg if (!has_array_of_type) { 457ff1e7120SSebastian Grimberg // Allocate if array is not yet allocated 458ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, mem_type, CEED_COPY_VALUES, NULL)); 459ff1e7120SSebastian Grimberg } else { 460ff1e7120SSebastian Grimberg // Select dirty array 461ff1e7120SSebastian Grimberg switch (mem_type) { 462ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 463ff1e7120SSebastian Grimberg if (impl->h_array_borrowed) impl->h_array = impl->h_array_borrowed; 464ff1e7120SSebastian Grimberg else impl->h_array = impl->h_array_owned; 465ff1e7120SSebastian Grimberg break; 466ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 467ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) impl->d_array = impl->d_array_borrowed; 468ff1e7120SSebastian Grimberg else impl->d_array = impl->d_array_owned; 469ff1e7120SSebastian Grimberg } 470ff1e7120SSebastian Grimberg } 471ff1e7120SSebastian Grimberg return CeedVectorGetArray_Cuda(vec, mem_type, array); 472ff1e7120SSebastian Grimberg } 473ff1e7120SSebastian Grimberg 474ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 475ff1e7120SSebastian Grimberg // Get the norm of a CeedVector 476ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 477ff1e7120SSebastian Grimberg static int CeedVectorNorm_Cuda(CeedVector vec, CeedNormType type, CeedScalar *norm) { 478ff1e7120SSebastian Grimberg Ceed ceed; 479ca735530SJeremy L Thompson CeedSize length; 480e84c3ebcSJeremy L Thompson #if (CUDA_VERSION < 12000) 481672b0f2aSSebastian Grimberg CeedSize num_calls; 482e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 483ca735530SJeremy L Thompson const CeedScalar *d_array; 484ca735530SJeremy L Thompson CeedVector_Cuda *impl; 485672b0f2aSSebastian Grimberg cublasHandle_t handle; 486ca735530SJeremy L Thompson 487ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 488ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 489ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 490ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetCublasHandle_Cuda(ceed, &handle)); 491ff1e7120SSebastian Grimberg 492e84c3ebcSJeremy L Thompson #if (CUDA_VERSION < 12000) 493f7c1b517Snbeams // With CUDA 12, we can use the 64-bit integer interface. Prior to that, 494f7c1b517Snbeams // we need to check if the vector is too long to handle with int32, 495b2165e7aSSebastian Grimberg // and if so, divide it into subsections for repeated cuBLAS calls. 496672b0f2aSSebastian Grimberg num_calls = length / INT_MAX; 497f7c1b517Snbeams if (length % INT_MAX > 0) num_calls += 1; 498e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 499f7c1b517Snbeams 500ff1e7120SSebastian Grimberg // Compute norm 501ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &d_array)); 502ff1e7120SSebastian Grimberg switch (type) { 503ff1e7120SSebastian Grimberg case CEED_NORM_1: { 504f6f49adbSnbeams *norm = 0.0; 505e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 506e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) // We have CUDA 12, and can use 64-bit integers 507f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 508e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 509f7c1b517Snbeams float sub_norm = 0.0; 510f7c1b517Snbeams float *d_array_start; 511ca735530SJeremy L Thompson 512f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 513f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 514f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 515f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 516ca735530SJeremy L Thompson 517f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 518f7c1b517Snbeams *norm += sub_norm; 519f7c1b517Snbeams } 520e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 521e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 522e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 523f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 524e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 525f7c1b517Snbeams double sub_norm = 0.0; 526f7c1b517Snbeams double *d_array_start; 527ca735530SJeremy L Thompson 528f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 529f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 530f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 531f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 532ca735530SJeremy L Thompson 533f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 534f7c1b517Snbeams *norm += sub_norm; 535f7c1b517Snbeams } 536e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 537e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 538ff1e7120SSebastian Grimberg break; 539ff1e7120SSebastian Grimberg } 540ff1e7120SSebastian Grimberg case CEED_NORM_2: { 541e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 542e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 543f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 544e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 545f7c1b517Snbeams float sub_norm = 0.0, norm_sum = 0.0; 546f7c1b517Snbeams float *d_array_start; 547ca735530SJeremy L Thompson 548f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 549f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 550f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 551f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 552ca735530SJeremy L Thompson 553f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 554f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 555f7c1b517Snbeams } 556f7c1b517Snbeams *norm = sqrt(norm_sum); 557e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 558e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 559e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 560f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 561e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 562f7c1b517Snbeams double sub_norm = 0.0, norm_sum = 0.0; 563f7c1b517Snbeams double *d_array_start; 564ca735530SJeremy L Thompson 565f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 566f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 567f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 568f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 569ca735530SJeremy L Thompson 570f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 571f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 572f7c1b517Snbeams } 573f7c1b517Snbeams *norm = sqrt(norm_sum); 574e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 575e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 576ff1e7120SSebastian Grimberg break; 577ff1e7120SSebastian Grimberg } 578ff1e7120SSebastian Grimberg case CEED_NORM_MAX: { 579e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 580e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 581ca735530SJeremy L Thompson int64_t index; 582ca735530SJeremy L Thompson CeedScalar norm_no_abs; 583ca735530SJeremy L Thompson 584ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax_64(handle, (int64_t)length, (float *)d_array, 1, &index)); 585ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 586ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 587e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 588ca735530SJeremy L Thompson CeedInt index; 589f7c1b517Snbeams float sub_max = 0.0, current_max = 0.0; 590f7c1b517Snbeams float *d_array_start; 591ca735530SJeremy L Thompson 592f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 593f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 594f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 595f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 596ca735530SJeremy L Thompson 597ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &index)); 598ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 599f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 600f7c1b517Snbeams } 601f7c1b517Snbeams *norm = current_max; 602e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 603e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 604e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 605ca735530SJeremy L Thompson int64_t index; 606ca735530SJeremy L Thompson CeedScalar norm_no_abs; 607ca735530SJeremy L Thompson 608ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax_64(handle, (int64_t)length, (double *)d_array, 1, &index)); 609ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 610ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 611e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 612ca735530SJeremy L Thompson CeedInt index; 613f7c1b517Snbeams double sub_max = 0.0, current_max = 0.0; 614f7c1b517Snbeams double *d_array_start; 615ca735530SJeremy L Thompson 616f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 617f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 618f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 619f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 620ca735530SJeremy L Thompson 621ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &index)); 622ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 623f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 624f7c1b517Snbeams } 625f7c1b517Snbeams *norm = current_max; 626e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 627e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 628ff1e7120SSebastian Grimberg break; 629ff1e7120SSebastian Grimberg } 630ff1e7120SSebastian Grimberg } 631ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(vec, &d_array)); 6329bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 633ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 634ff1e7120SSebastian Grimberg } 635ff1e7120SSebastian Grimberg 636ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 637ff1e7120SSebastian Grimberg // Take reciprocal of a vector on host 638ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 639f7c1b517Snbeams static int CeedHostReciprocal_Cuda(CeedScalar *h_array, CeedSize length) { 640f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) { 641ff1e7120SSebastian Grimberg if (fabs(h_array[i]) > CEED_EPSILON) h_array[i] = 1. / h_array[i]; 642ff1e7120SSebastian Grimberg } 643ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 644ff1e7120SSebastian Grimberg } 645ff1e7120SSebastian Grimberg 646ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 647ff1e7120SSebastian Grimberg // Take reciprocal of a vector on device (impl in .cu file) 648ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 649f7c1b517Snbeams int CeedDeviceReciprocal_Cuda(CeedScalar *d_array, CeedSize length); 650ff1e7120SSebastian Grimberg 651ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 652ff1e7120SSebastian Grimberg // Take reciprocal of a vector 653ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 654ff1e7120SSebastian Grimberg static int CeedVectorReciprocal_Cuda(CeedVector vec) { 655ff1e7120SSebastian Grimberg CeedSize length; 656ca735530SJeremy L Thompson CeedVector_Cuda *impl; 657ff1e7120SSebastian Grimberg 658ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 659ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 660ff1e7120SSebastian Grimberg // Set value for synced device/host array 661ff1e7120SSebastian Grimberg if (impl->d_array) CeedCallBackend(CeedDeviceReciprocal_Cuda(impl->d_array, length)); 662ff1e7120SSebastian Grimberg if (impl->h_array) CeedCallBackend(CeedHostReciprocal_Cuda(impl->h_array, length)); 663ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 664ff1e7120SSebastian Grimberg } 665ff1e7120SSebastian Grimberg 666ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 667ff1e7120SSebastian Grimberg // Compute x = alpha x on the host 668ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 669f7c1b517Snbeams static int CeedHostScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length) { 670f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) x_array[i] *= alpha; 671ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 672ff1e7120SSebastian Grimberg } 673ff1e7120SSebastian Grimberg 674ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 675ff1e7120SSebastian Grimberg // Compute x = alpha x on device (impl in .cu file) 676ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 677f7c1b517Snbeams int CeedDeviceScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length); 678ff1e7120SSebastian Grimberg 679ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 680ff1e7120SSebastian Grimberg // Compute x = alpha x 681ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 682ff1e7120SSebastian Grimberg static int CeedVectorScale_Cuda(CeedVector x, CeedScalar alpha) { 683ff1e7120SSebastian Grimberg CeedSize length; 684e84c3ebcSJeremy L Thompson CeedVector_Cuda *impl; 685ff1e7120SSebastian Grimberg 686e84c3ebcSJeremy L Thompson CeedCallBackend(CeedVectorGetData(x, &impl)); 687ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(x, &length)); 688ff1e7120SSebastian Grimberg // Set value for synced device/host array 689e84c3ebcSJeremy L Thompson if (impl->d_array) { 690e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 691e84c3ebcSJeremy L Thompson cublasHandle_t handle; 692e84c3ebcSJeremy L Thompson 693e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(CeedVectorReturnCeed(x), &handle)); 694e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 695e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(x), cublasSscal_64(handle, (int64_t)length, &alpha, impl->d_array, 1)); 696e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 697e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(x), cublasDscal_64(handle, (int64_t)length, &alpha, impl->d_array, 1)); 698e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 699e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 700e84c3ebcSJeremy L Thompson CeedCallBackend(CeedDeviceScale_Cuda(impl->d_array, alpha, length)); 701e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 702e84c3ebcSJeremy L Thompson impl->h_array = NULL; 703e84c3ebcSJeremy L Thompson } else if (impl->h_array) { 704e84c3ebcSJeremy L Thompson CeedCallBackend(CeedHostScale_Cuda(impl->h_array, alpha, length)); 705e84c3ebcSJeremy L Thompson impl->d_array = NULL; 706e84c3ebcSJeremy L Thompson } 707ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 708ff1e7120SSebastian Grimberg } 709ff1e7120SSebastian Grimberg 710ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 711ff1e7120SSebastian Grimberg // Compute y = alpha x + y on the host 712ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 713f7c1b517Snbeams static int CeedHostAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length) { 714f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) y_array[i] += alpha * x_array[i]; 715ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 716ff1e7120SSebastian Grimberg } 717ff1e7120SSebastian Grimberg 718ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 719ff1e7120SSebastian Grimberg // Compute y = alpha x + y on device (impl in .cu file) 720ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 721f7c1b517Snbeams int CeedDeviceAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length); 722ff1e7120SSebastian Grimberg 723ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 724ff1e7120SSebastian Grimberg // Compute y = alpha x + y 725ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 726ff1e7120SSebastian Grimberg static int CeedVectorAXPY_Cuda(CeedVector y, CeedScalar alpha, CeedVector x) { 727ca735530SJeremy L Thompson CeedSize length; 728ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 729ca735530SJeremy L Thompson 730ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 731ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 732ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 733ff1e7120SSebastian Grimberg // Set value for synced device/host array 734ff1e7120SSebastian Grimberg if (y_impl->d_array) { 735ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 736e84c3ebcSJeremy L Thompson #if (CUDA_VERSION >= 12000) 737e84c3ebcSJeremy L Thompson cublasHandle_t handle; 738e84c3ebcSJeremy L Thompson 739e84c3ebcSJeremy L Thompson CeedCallBackend(CeedGetCublasHandle_Cuda(CeedVectorReturnCeed(y), &handle)); 740e84c3ebcSJeremy L Thompson #if defined(CEED_SCALAR_IS_FP32) 741e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(y), cublasSaxpy_64(handle, (int64_t)length, &alpha, x_impl->d_array, 1, y_impl->d_array, 1)); 742e84c3ebcSJeremy L Thompson #else /* CEED_SCALAR */ 743e84c3ebcSJeremy L Thompson CeedCallCublas(CeedVectorReturnCeed(y), cublasDaxpy_64(handle, (int64_t)length, &alpha, x_impl->d_array, 1, y_impl->d_array, 1)); 744e84c3ebcSJeremy L Thompson #endif /* CEED_SCALAR */ 745e84c3ebcSJeremy L Thompson #else /* CUDA_VERSION */ 746ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPY_Cuda(y_impl->d_array, alpha, x_impl->d_array, length)); 747e84c3ebcSJeremy L Thompson #endif /* CUDA_VERSION */ 748e84c3ebcSJeremy L Thompson y_impl->h_array = NULL; 749e84c3ebcSJeremy L Thompson } else if (y_impl->h_array) { 750ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 751ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPY_Cuda(y_impl->h_array, alpha, x_impl->h_array, length)); 752e84c3ebcSJeremy L Thompson y_impl->d_array = NULL; 753ff1e7120SSebastian Grimberg } 754ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 755ff1e7120SSebastian Grimberg } 756ff1e7120SSebastian Grimberg 757ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 758ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on the host 759ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 760f7c1b517Snbeams static int CeedHostAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length) { 761aa67b842SZach Atkins for (CeedSize i = 0; i < length; i++) y_array[i] = alpha * x_array[i] + beta * y_array[i]; 762ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 763ff1e7120SSebastian Grimberg } 764ff1e7120SSebastian Grimberg 765ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 766ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on device (impl in .cu file) 767ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 768f7c1b517Snbeams int CeedDeviceAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length); 769ff1e7120SSebastian Grimberg 770ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 771ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y 772ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 773ff1e7120SSebastian Grimberg static int CeedVectorAXPBY_Cuda(CeedVector y, CeedScalar alpha, CeedScalar beta, CeedVector x) { 774ca735530SJeremy L Thompson CeedSize length; 775ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 776ca735530SJeremy L Thompson 777ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 778ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 779ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 780ff1e7120SSebastian Grimberg // Set value for synced device/host array 781ff1e7120SSebastian Grimberg if (y_impl->d_array) { 782ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 783ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPBY_Cuda(y_impl->d_array, alpha, beta, x_impl->d_array, length)); 784ff1e7120SSebastian Grimberg } 785ff1e7120SSebastian Grimberg if (y_impl->h_array) { 786ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 787ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPBY_Cuda(y_impl->h_array, alpha, beta, x_impl->h_array, length)); 788ff1e7120SSebastian Grimberg } 789ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 790ff1e7120SSebastian Grimberg } 791ff1e7120SSebastian Grimberg 792ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 793ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on the host 794ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 795f7c1b517Snbeams static int CeedHostPointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length) { 796f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) w_array[i] = x_array[i] * y_array[i]; 797ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 798ff1e7120SSebastian Grimberg } 799ff1e7120SSebastian Grimberg 800ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 801ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on device (impl in .cu file) 802ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 803f7c1b517Snbeams int CeedDevicePointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length); 804ff1e7120SSebastian Grimberg 805ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 806ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y 807ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 808ff1e7120SSebastian Grimberg static int CeedVectorPointwiseMult_Cuda(CeedVector w, CeedVector x, CeedVector y) { 809ca735530SJeremy L Thompson CeedSize length; 810ff1e7120SSebastian Grimberg CeedVector_Cuda *w_impl, *x_impl, *y_impl; 811ca735530SJeremy L Thompson 812ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(w, &w_impl)); 813ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 814ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 815ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(w, &length)); 816ff1e7120SSebastian Grimberg // Set value for synced device/host array 817ff1e7120SSebastian Grimberg if (!w_impl->d_array && !w_impl->h_array) { 818ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetValue(w, 0.0)); 819ff1e7120SSebastian Grimberg } 820ff1e7120SSebastian Grimberg if (w_impl->d_array) { 821ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 822ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_DEVICE)); 823ff1e7120SSebastian Grimberg CeedCallBackend(CeedDevicePointwiseMult_Cuda(w_impl->d_array, x_impl->d_array, y_impl->d_array, length)); 824ff1e7120SSebastian Grimberg } 825ff1e7120SSebastian Grimberg if (w_impl->h_array) { 826ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 827ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_HOST)); 828ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostPointwiseMult_Cuda(w_impl->h_array, x_impl->h_array, y_impl->h_array, length)); 829ff1e7120SSebastian Grimberg } 830ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 831ff1e7120SSebastian Grimberg } 832ff1e7120SSebastian Grimberg 833ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 834ff1e7120SSebastian Grimberg // Destroy the vector 835ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 836ff1e7120SSebastian Grimberg static int CeedVectorDestroy_Cuda(const CeedVector vec) { 837ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 838ff1e7120SSebastian Grimberg 839ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 8406e536b99SJeremy L Thompson CeedCallCuda(CeedVectorReturnCeed(vec), cudaFree(impl->d_array_owned)); 841ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_array_owned)); 842ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 843ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 844ff1e7120SSebastian Grimberg } 845ff1e7120SSebastian Grimberg 846ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 847ff1e7120SSebastian Grimberg // Create a vector of the specified length (does not allocate memory) 848ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 849ff1e7120SSebastian Grimberg int CeedVectorCreate_Cuda(CeedSize n, CeedVector vec) { 850ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 851ff1e7120SSebastian Grimberg Ceed ceed; 852ff1e7120SSebastian Grimberg 853ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 854ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasValidArray", CeedVectorHasValidArray_Cuda)); 855ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasBorrowedArrayOfType", CeedVectorHasBorrowedArrayOfType_Cuda)); 856ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetArray", CeedVectorSetArray_Cuda)); 857ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "TakeArray", CeedVectorTakeArray_Cuda)); 8583e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "CopyStrided", CeedVectorCopyStrided_Cuda)); 8593e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetValue", CeedVectorSetValue_Cuda)); 8603e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetValueStrided", CeedVectorSetValueStrided_Cuda)); 861ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SyncArray", CeedVectorSyncArray_Cuda)); 862ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArray", CeedVectorGetArray_Cuda)); 863ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayRead", CeedVectorGetArrayRead_Cuda)); 864ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayWrite", CeedVectorGetArrayWrite_Cuda)); 865ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Norm", CeedVectorNorm_Cuda)); 866ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Reciprocal", CeedVectorReciprocal_Cuda)); 8673e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Scale", CeedVectorScale_Cuda)); 8683e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPY", CeedVectorAXPY_Cuda)); 8693e961e14SJeremy L Thompson CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPBY", CeedVectorAXPBY_Cuda)); 870ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "PointwiseMult", CeedVectorPointwiseMult_Cuda)); 871ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Destroy", CeedVectorDestroy_Cuda)); 8729bc66399SJeremy L Thompson CeedCallBackend(CeedDestroy(&ceed)); 873ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 874ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetData(vec, impl)); 875ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 876ff1e7120SSebastian Grimberg } 877ff1e7120SSebastian Grimberg 878ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 879