1ff1e7120SSebastian Grimberg // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2ff1e7120SSebastian Grimberg // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3ff1e7120SSebastian Grimberg // 4ff1e7120SSebastian Grimberg // SPDX-License-Identifier: BSD-2-Clause 5ff1e7120SSebastian Grimberg // 6ff1e7120SSebastian Grimberg // This file is part of CEED: http://github.com/ceed 7ff1e7120SSebastian Grimberg 8ff1e7120SSebastian Grimberg #include <ceed.h> 9ff1e7120SSebastian Grimberg #include <ceed/backend.h> 10ff1e7120SSebastian Grimberg #include <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) { 42ff1e7120SSebastian Grimberg Ceed ceed; 43ca735530SJeremy L Thompson CeedSize length; 44672b0f2aSSebastian Grimberg size_t bytes; 45ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 46ca735530SJeremy L Thompson 47ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 48ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 49ff1e7120SSebastian Grimberg 50ff1e7120SSebastian Grimberg CeedCheck(impl->h_array, ceed, CEED_ERROR_BACKEND, "No valid host data to sync to device"); 51ff1e7120SSebastian Grimberg 52ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 53672b0f2aSSebastian Grimberg bytes = length * sizeof(CeedScalar); 54ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) { 55ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_borrowed; 56ff1e7120SSebastian Grimberg } else if (impl->d_array_owned) { 57ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 58ff1e7120SSebastian Grimberg } else { 59ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_array_owned, bytes)); 60ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 61ff1e7120SSebastian Grimberg } 62ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMemcpy(impl->d_array, impl->h_array, bytes, cudaMemcpyHostToDevice)); 63ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 64ff1e7120SSebastian Grimberg } 65ff1e7120SSebastian Grimberg 66ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 67ff1e7120SSebastian Grimberg // Sync device to host 68ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 69ff1e7120SSebastian Grimberg static inline int CeedVectorSyncD2H_Cuda(const CeedVector vec) { 70ff1e7120SSebastian Grimberg Ceed ceed; 71ca735530SJeremy L Thompson CeedSize length; 72ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 73ca735530SJeremy L Thompson 74ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 75ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 76ff1e7120SSebastian Grimberg 77ff1e7120SSebastian Grimberg CeedCheck(impl->d_array, ceed, CEED_ERROR_BACKEND, "No valid device data to sync to host"); 78ff1e7120SSebastian Grimberg 79ff1e7120SSebastian Grimberg if (impl->h_array_borrowed) { 80ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_borrowed; 81ff1e7120SSebastian Grimberg } else if (impl->h_array_owned) { 82ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 83ff1e7120SSebastian Grimberg } else { 84ff1e7120SSebastian Grimberg CeedSize length; 85ca735530SJeremy L Thompson 86ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 87ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(length, &impl->h_array_owned)); 88ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 89ff1e7120SSebastian Grimberg } 90ff1e7120SSebastian Grimberg 91ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 92ff1e7120SSebastian Grimberg size_t bytes = length * sizeof(CeedScalar); 93ff1e7120SSebastian Grimberg 94ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(impl->h_array, impl->d_array, bytes, cudaMemcpyDeviceToHost)); 95ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 96ff1e7120SSebastian Grimberg } 97ff1e7120SSebastian Grimberg 98ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 99ff1e7120SSebastian Grimberg // Sync arrays 100ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 101ff1e7120SSebastian Grimberg static int CeedVectorSyncArray_Cuda(const CeedVector vec, CeedMemType mem_type) { 102ff1e7120SSebastian Grimberg bool need_sync = false; 103ca735530SJeremy L Thompson 104ca735530SJeremy L Thompson // Check whether device/host sync is needed 105ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorNeedSync_Cuda(vec, mem_type, &need_sync)); 106ff1e7120SSebastian Grimberg if (!need_sync) return CEED_ERROR_SUCCESS; 107ff1e7120SSebastian Grimberg 108ff1e7120SSebastian Grimberg switch (mem_type) { 109ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 110ff1e7120SSebastian Grimberg return CeedVectorSyncD2H_Cuda(vec); 111ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 112ff1e7120SSebastian Grimberg return CeedVectorSyncH2D_Cuda(vec); 113ff1e7120SSebastian Grimberg } 114ff1e7120SSebastian Grimberg return CEED_ERROR_UNSUPPORTED; 115ff1e7120SSebastian Grimberg } 116ff1e7120SSebastian Grimberg 117ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 118ff1e7120SSebastian Grimberg // Set all pointers as invalid 119ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 120ff1e7120SSebastian Grimberg static inline int CeedVectorSetAllInvalid_Cuda(const CeedVector vec) { 121ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 122ff1e7120SSebastian Grimberg 123ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 124ff1e7120SSebastian Grimberg impl->h_array = NULL; 125ff1e7120SSebastian Grimberg impl->d_array = NULL; 126ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 127ff1e7120SSebastian Grimberg } 128ff1e7120SSebastian Grimberg 129ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 130ff1e7120SSebastian Grimberg // Check if CeedVector has any valid pointer 131ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 132ff1e7120SSebastian Grimberg static inline int CeedVectorHasValidArray_Cuda(const CeedVector vec, bool *has_valid_array) { 133ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 134ca735530SJeremy L Thompson 135ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 1361c66c397SJeremy L Thompson *has_valid_array = impl->h_array || impl->d_array; 137ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 138ff1e7120SSebastian Grimberg } 139ff1e7120SSebastian Grimberg 140ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 141ff1e7120SSebastian Grimberg // Check if has array of given type 142ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 143ff1e7120SSebastian Grimberg static inline int CeedVectorHasArrayOfType_Cuda(const CeedVector vec, CeedMemType mem_type, bool *has_array_of_type) { 144ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 145ff1e7120SSebastian Grimberg 146ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 147ff1e7120SSebastian Grimberg switch (mem_type) { 148ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 1491c66c397SJeremy L Thompson *has_array_of_type = impl->h_array_borrowed || impl->h_array_owned; 150ff1e7120SSebastian Grimberg break; 151ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 1521c66c397SJeremy L Thompson *has_array_of_type = impl->d_array_borrowed || impl->d_array_owned; 153ff1e7120SSebastian Grimberg break; 154ff1e7120SSebastian Grimberg } 155ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 156ff1e7120SSebastian Grimberg } 157ff1e7120SSebastian Grimberg 158ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 159ff1e7120SSebastian Grimberg // Check if has borrowed array of given type 160ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 161ff1e7120SSebastian Grimberg static inline int CeedVectorHasBorrowedArrayOfType_Cuda(const CeedVector vec, CeedMemType mem_type, bool *has_borrowed_array_of_type) { 162ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 163ff1e7120SSebastian Grimberg 164ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 165ff1e7120SSebastian Grimberg switch (mem_type) { 166ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 1671c66c397SJeremy L Thompson *has_borrowed_array_of_type = impl->h_array_borrowed; 168ff1e7120SSebastian Grimberg break; 169ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 1701c66c397SJeremy L Thompson *has_borrowed_array_of_type = impl->d_array_borrowed; 171ff1e7120SSebastian Grimberg break; 172ff1e7120SSebastian Grimberg } 173ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 174ff1e7120SSebastian Grimberg } 175ff1e7120SSebastian Grimberg 176ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 177ff1e7120SSebastian Grimberg // Set array from host 178ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 179ff1e7120SSebastian Grimberg static int CeedVectorSetArrayHost_Cuda(const CeedVector vec, const CeedCopyMode copy_mode, CeedScalar *array) { 180ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 181ff1e7120SSebastian Grimberg 182ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 183ff1e7120SSebastian Grimberg switch (copy_mode) { 184ff1e7120SSebastian Grimberg case CEED_COPY_VALUES: { 185ff1e7120SSebastian Grimberg if (!impl->h_array_owned) { 186ca735530SJeremy L Thompson CeedSize length; 187ca735530SJeremy L Thompson 188ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 189ff1e7120SSebastian Grimberg CeedCallBackend(CeedMalloc(length, &impl->h_array_owned)); 190ff1e7120SSebastian Grimberg } 191ff1e7120SSebastian Grimberg impl->h_array_borrowed = NULL; 192ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 193ff1e7120SSebastian Grimberg if (array) { 194ff1e7120SSebastian Grimberg CeedSize length; 195672b0f2aSSebastian Grimberg size_t bytes; 196ca735530SJeremy L Thompson 197ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 198672b0f2aSSebastian Grimberg bytes = length * sizeof(CeedScalar); 199ff1e7120SSebastian Grimberg memcpy(impl->h_array, array, bytes); 200ff1e7120SSebastian Grimberg } 201ff1e7120SSebastian Grimberg } break; 202ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 203ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_array_owned)); 204ff1e7120SSebastian Grimberg impl->h_array_owned = array; 205ff1e7120SSebastian Grimberg impl->h_array_borrowed = NULL; 206ff1e7120SSebastian Grimberg impl->h_array = array; 207ff1e7120SSebastian Grimberg break; 208ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 209ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_array_owned)); 210ff1e7120SSebastian Grimberg impl->h_array_borrowed = array; 211ff1e7120SSebastian Grimberg impl->h_array = array; 212ff1e7120SSebastian Grimberg break; 213ff1e7120SSebastian Grimberg } 214ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 215ff1e7120SSebastian Grimberg } 216ff1e7120SSebastian Grimberg 217ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 218ff1e7120SSebastian Grimberg // Set array from device 219ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 220ff1e7120SSebastian Grimberg static int CeedVectorSetArrayDevice_Cuda(const CeedVector vec, const CeedCopyMode copy_mode, CeedScalar *array) { 221ff1e7120SSebastian Grimberg Ceed ceed; 222ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 223ff1e7120SSebastian Grimberg 224ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 225ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 226ff1e7120SSebastian Grimberg switch (copy_mode) { 227ff1e7120SSebastian Grimberg case CEED_COPY_VALUES: { 228ff1e7120SSebastian Grimberg CeedSize length; 229672b0f2aSSebastian Grimberg size_t bytes; 230ca735530SJeremy L Thompson 231ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(vec, &length)); 232672b0f2aSSebastian Grimberg bytes = length * sizeof(CeedScalar); 233ff1e7120SSebastian Grimberg if (!impl->d_array_owned) { 234ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaMalloc((void **)&impl->d_array_owned, bytes)); 235ff1e7120SSebastian Grimberg } 236b2165e7aSSebastian Grimberg impl->d_array_borrowed = NULL; 2373ce2313bSJeremy L Thompson impl->d_array = impl->d_array_owned; 238ff1e7120SSebastian Grimberg if (array) CeedCallCuda(ceed, cudaMemcpy(impl->d_array, array, bytes, cudaMemcpyDeviceToDevice)); 239ff1e7120SSebastian Grimberg } break; 240ff1e7120SSebastian Grimberg case CEED_OWN_POINTER: 241ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_array_owned)); 242ff1e7120SSebastian Grimberg impl->d_array_owned = array; 243ff1e7120SSebastian Grimberg impl->d_array_borrowed = NULL; 244ff1e7120SSebastian Grimberg impl->d_array = array; 245ff1e7120SSebastian Grimberg break; 246ff1e7120SSebastian Grimberg case CEED_USE_POINTER: 247ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_array_owned)); 248ff1e7120SSebastian Grimberg impl->d_array_owned = NULL; 249ff1e7120SSebastian Grimberg impl->d_array_borrowed = array; 250ff1e7120SSebastian Grimberg impl->d_array = array; 251ff1e7120SSebastian Grimberg break; 252ff1e7120SSebastian Grimberg } 253ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 254ff1e7120SSebastian Grimberg } 255ff1e7120SSebastian Grimberg 256ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 257ff1e7120SSebastian Grimberg // Set the array used by a vector, 258ff1e7120SSebastian Grimberg // freeing any previously allocated array if applicable 259ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 260ff1e7120SSebastian Grimberg static int CeedVectorSetArray_Cuda(const CeedVector vec, const CeedMemType mem_type, const CeedCopyMode copy_mode, CeedScalar *array) { 261ff1e7120SSebastian Grimberg Ceed ceed; 262ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 263ff1e7120SSebastian Grimberg 264ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 265ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 266ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetAllInvalid_Cuda(vec)); 267ff1e7120SSebastian Grimberg switch (mem_type) { 268ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 269ff1e7120SSebastian Grimberg return CeedVectorSetArrayHost_Cuda(vec, copy_mode, array); 270ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 271ff1e7120SSebastian Grimberg return CeedVectorSetArrayDevice_Cuda(vec, copy_mode, array); 272ff1e7120SSebastian Grimberg } 273ff1e7120SSebastian Grimberg return CEED_ERROR_UNSUPPORTED; 274ff1e7120SSebastian Grimberg } 275ff1e7120SSebastian Grimberg 276ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 277ff1e7120SSebastian Grimberg // Set host array to value 278ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 279f7c1b517Snbeams static int CeedHostSetValue_Cuda(CeedScalar *h_array, CeedSize length, CeedScalar val) { 280f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) h_array[i] = val; 281ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 282ff1e7120SSebastian Grimberg } 283ff1e7120SSebastian Grimberg 284ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 285ff1e7120SSebastian Grimberg // Set device array to value (impl in .cu file) 286ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 287f7c1b517Snbeams int CeedDeviceSetValue_Cuda(CeedScalar *d_array, CeedSize length, CeedScalar val); 288ff1e7120SSebastian Grimberg 289ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 290f7c1b517Snbeams // Set a vector to a value 291ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 292ff1e7120SSebastian Grimberg static int CeedVectorSetValue_Cuda(CeedVector vec, CeedScalar val) { 293ff1e7120SSebastian Grimberg Ceed ceed; 294ff1e7120SSebastian Grimberg CeedSize length; 295ca735530SJeremy L Thompson CeedVector_Cuda *impl; 296ff1e7120SSebastian Grimberg 297ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 298ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 299ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 300ff1e7120SSebastian Grimberg // Set value for synced device/host array 301ff1e7120SSebastian Grimberg if (!impl->d_array && !impl->h_array) { 302ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) { 303ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_borrowed; 304ff1e7120SSebastian Grimberg } else if (impl->h_array_borrowed) { 305ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_borrowed; 306ff1e7120SSebastian Grimberg } else if (impl->d_array_owned) { 307ff1e7120SSebastian Grimberg impl->d_array = impl->d_array_owned; 308ff1e7120SSebastian Grimberg } else if (impl->h_array_owned) { 309ff1e7120SSebastian Grimberg impl->h_array = impl->h_array_owned; 310ff1e7120SSebastian Grimberg } else { 311ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, CEED_MEM_DEVICE, CEED_COPY_VALUES, NULL)); 312ff1e7120SSebastian Grimberg } 313ff1e7120SSebastian Grimberg } 314ff1e7120SSebastian Grimberg if (impl->d_array) { 315ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceSetValue_Cuda(impl->d_array, length, val)); 316ff1e7120SSebastian Grimberg impl->h_array = NULL; 317ff1e7120SSebastian Grimberg } 318ff1e7120SSebastian Grimberg if (impl->h_array) { 319ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostSetValue_Cuda(impl->h_array, length, val)); 320ff1e7120SSebastian Grimberg impl->d_array = NULL; 321ff1e7120SSebastian Grimberg } 322ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 323ff1e7120SSebastian Grimberg } 324ff1e7120SSebastian Grimberg 325ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 326ff1e7120SSebastian Grimberg // Vector Take Array 327ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 328ff1e7120SSebastian Grimberg static int CeedVectorTakeArray_Cuda(CeedVector vec, CeedMemType mem_type, CeedScalar **array) { 329ff1e7120SSebastian Grimberg Ceed ceed; 330ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 331ff1e7120SSebastian Grimberg 332ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 333ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 334ff1e7120SSebastian Grimberg // Sync array to requested mem_type 335ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 336ff1e7120SSebastian Grimberg // Update pointer 337ff1e7120SSebastian Grimberg switch (mem_type) { 338ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 339ff1e7120SSebastian Grimberg (*array) = impl->h_array_borrowed; 340ff1e7120SSebastian Grimberg impl->h_array_borrowed = NULL; 341ff1e7120SSebastian Grimberg impl->h_array = NULL; 342ff1e7120SSebastian Grimberg break; 343ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 344ff1e7120SSebastian Grimberg (*array) = impl->d_array_borrowed; 345ff1e7120SSebastian Grimberg impl->d_array_borrowed = NULL; 346ff1e7120SSebastian Grimberg impl->d_array = NULL; 347ff1e7120SSebastian Grimberg break; 348ff1e7120SSebastian Grimberg } 349ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 350ff1e7120SSebastian Grimberg } 351ff1e7120SSebastian Grimberg 352ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 353ff1e7120SSebastian Grimberg // Core logic for array syncronization for GetArray. 354ff1e7120SSebastian Grimberg // If a different memory type is most up to date, this will perform a copy 355ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 356ff1e7120SSebastian Grimberg static int CeedVectorGetArrayCore_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 357ff1e7120SSebastian Grimberg Ceed ceed; 358ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 359ff1e7120SSebastian Grimberg 360ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 361ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 362ff1e7120SSebastian Grimberg // Sync array to requested mem_type 363ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(vec, mem_type)); 364ff1e7120SSebastian Grimberg // Update pointer 365ff1e7120SSebastian Grimberg switch (mem_type) { 366ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 367ff1e7120SSebastian Grimberg *array = impl->h_array; 368ff1e7120SSebastian Grimberg break; 369ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 370ff1e7120SSebastian Grimberg *array = impl->d_array; 371ff1e7120SSebastian Grimberg break; 372ff1e7120SSebastian Grimberg } 373ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 374ff1e7120SSebastian Grimberg } 375ff1e7120SSebastian Grimberg 376ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 377ff1e7120SSebastian Grimberg // Get read-only access to a vector via the specified mem_type 378ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 379ff1e7120SSebastian Grimberg static int CeedVectorGetArrayRead_Cuda(const CeedVector vec, const CeedMemType mem_type, const CeedScalar **array) { 380ff1e7120SSebastian Grimberg return CeedVectorGetArrayCore_Cuda(vec, mem_type, (CeedScalar **)array); 381ff1e7120SSebastian Grimberg } 382ff1e7120SSebastian Grimberg 383ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 384ff1e7120SSebastian Grimberg // Get read/write access to a vector via the specified mem_type 385ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 386ff1e7120SSebastian Grimberg static int CeedVectorGetArray_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 387ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 388ca735530SJeremy L Thompson 389ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(vec, &impl)); 390ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayCore_Cuda(vec, mem_type, array)); 391ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetAllInvalid_Cuda(vec)); 392ff1e7120SSebastian Grimberg switch (mem_type) { 393ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 394ff1e7120SSebastian Grimberg impl->h_array = *array; 395ff1e7120SSebastian Grimberg break; 396ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 397ff1e7120SSebastian Grimberg impl->d_array = *array; 398ff1e7120SSebastian Grimberg break; 399ff1e7120SSebastian Grimberg } 400ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 401ff1e7120SSebastian Grimberg } 402ff1e7120SSebastian Grimberg 403ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 404ff1e7120SSebastian Grimberg // Get write access to a vector via the specified mem_type 405ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 406ff1e7120SSebastian Grimberg static int CeedVectorGetArrayWrite_Cuda(const CeedVector vec, const CeedMemType mem_type, CeedScalar **array) { 407ff1e7120SSebastian Grimberg bool has_array_of_type = true; 408ca735530SJeremy L Thompson CeedVector_Cuda *impl; 409ca735530SJeremy L Thompson 410ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 411ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorHasArrayOfType_Cuda(vec, mem_type, &has_array_of_type)); 412ff1e7120SSebastian Grimberg if (!has_array_of_type) { 413ff1e7120SSebastian Grimberg // Allocate if array is not yet allocated 414ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetArray(vec, mem_type, CEED_COPY_VALUES, NULL)); 415ff1e7120SSebastian Grimberg } else { 416ff1e7120SSebastian Grimberg // Select dirty array 417ff1e7120SSebastian Grimberg switch (mem_type) { 418ff1e7120SSebastian Grimberg case CEED_MEM_HOST: 419ff1e7120SSebastian Grimberg if (impl->h_array_borrowed) impl->h_array = impl->h_array_borrowed; 420ff1e7120SSebastian Grimberg else impl->h_array = impl->h_array_owned; 421ff1e7120SSebastian Grimberg break; 422ff1e7120SSebastian Grimberg case CEED_MEM_DEVICE: 423ff1e7120SSebastian Grimberg if (impl->d_array_borrowed) impl->d_array = impl->d_array_borrowed; 424ff1e7120SSebastian Grimberg else impl->d_array = impl->d_array_owned; 425ff1e7120SSebastian Grimberg } 426ff1e7120SSebastian Grimberg } 427ff1e7120SSebastian Grimberg return CeedVectorGetArray_Cuda(vec, mem_type, array); 428ff1e7120SSebastian Grimberg } 429ff1e7120SSebastian Grimberg 430ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 431ff1e7120SSebastian Grimberg // Get the norm of a CeedVector 432ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 433ff1e7120SSebastian Grimberg static int CeedVectorNorm_Cuda(CeedVector vec, CeedNormType type, CeedScalar *norm) { 434ff1e7120SSebastian Grimberg Ceed ceed; 435ca735530SJeremy L Thompson CeedSize length; 436672b0f2aSSebastian Grimberg #if CUDA_VERSION < 12000 437672b0f2aSSebastian Grimberg CeedSize num_calls; 438672b0f2aSSebastian Grimberg #endif 439ca735530SJeremy L Thompson const CeedScalar *d_array; 440ca735530SJeremy L Thompson CeedVector_Cuda *impl; 441672b0f2aSSebastian Grimberg cublasHandle_t handle; 442ca735530SJeremy L Thompson 443ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 444ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 445ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 446ff1e7120SSebastian Grimberg CeedCallBackend(CeedGetCublasHandle_Cuda(ceed, &handle)); 447ff1e7120SSebastian Grimberg 448f7c1b517Snbeams #if CUDA_VERSION < 12000 449f7c1b517Snbeams // With CUDA 12, we can use the 64-bit integer interface. Prior to that, 450f7c1b517Snbeams // we need to check if the vector is too long to handle with int32, 451b2165e7aSSebastian Grimberg // and if so, divide it into subsections for repeated cuBLAS calls. 452672b0f2aSSebastian Grimberg num_calls = length / INT_MAX; 453f7c1b517Snbeams if (length % INT_MAX > 0) num_calls += 1; 454f7c1b517Snbeams #endif 455f7c1b517Snbeams 456ff1e7120SSebastian Grimberg // Compute norm 457ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, &d_array)); 458ff1e7120SSebastian Grimberg switch (type) { 459ff1e7120SSebastian Grimberg case CEED_NORM_1: { 460f6f49adbSnbeams *norm = 0.0; 461ff1e7120SSebastian Grimberg if (CEED_SCALAR_TYPE == CEED_SCALAR_FP32) { 462f7c1b517Snbeams #if CUDA_VERSION >= 12000 // We have CUDA 12, and can use 64-bit integers 463f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 464f7c1b517Snbeams #else 465f7c1b517Snbeams float sub_norm = 0.0; 466f7c1b517Snbeams float *d_array_start; 467ca735530SJeremy L Thompson 468f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 469f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 470f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 471f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 472ca735530SJeremy L Thompson 473f7c1b517Snbeams CeedCallCublas(ceed, cublasSasum(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 474f7c1b517Snbeams *norm += sub_norm; 475f7c1b517Snbeams } 476f7c1b517Snbeams #endif 477ff1e7120SSebastian Grimberg } else { 478f7c1b517Snbeams #if CUDA_VERSION >= 12000 479f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 480f7c1b517Snbeams #else 481f7c1b517Snbeams double sub_norm = 0.0; 482f7c1b517Snbeams double *d_array_start; 483ca735530SJeremy L Thompson 484f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 485f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 486f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 487f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 488ca735530SJeremy L Thompson 489f7c1b517Snbeams CeedCallCublas(ceed, cublasDasum(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 490f7c1b517Snbeams *norm += sub_norm; 491f7c1b517Snbeams } 492f7c1b517Snbeams #endif 493ff1e7120SSebastian Grimberg } 494ff1e7120SSebastian Grimberg break; 495ff1e7120SSebastian Grimberg } 496ff1e7120SSebastian Grimberg case CEED_NORM_2: { 497ff1e7120SSebastian Grimberg if (CEED_SCALAR_TYPE == CEED_SCALAR_FP32) { 498f7c1b517Snbeams #if CUDA_VERSION >= 12000 499f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2_64(handle, (int64_t)length, (float *)d_array, 1, (float *)norm)); 500f7c1b517Snbeams #else 501f7c1b517Snbeams float sub_norm = 0.0, norm_sum = 0.0; 502f7c1b517Snbeams float *d_array_start; 503ca735530SJeremy L Thompson 504f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 505f7c1b517Snbeams d_array_start = (float *)d_array + (CeedSize)(i)*INT_MAX; 506f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 507f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 508ca735530SJeremy L Thompson 509f7c1b517Snbeams CeedCallCublas(ceed, cublasSnrm2(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &sub_norm)); 510f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 511f7c1b517Snbeams } 512f7c1b517Snbeams *norm = sqrt(norm_sum); 513f7c1b517Snbeams #endif 514ff1e7120SSebastian Grimberg } else { 515f7c1b517Snbeams #if CUDA_VERSION >= 12000 516f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2_64(handle, (int64_t)length, (double *)d_array, 1, (double *)norm)); 517f7c1b517Snbeams #else 518f7c1b517Snbeams double sub_norm = 0.0, norm_sum = 0.0; 519f7c1b517Snbeams double *d_array_start; 520ca735530SJeremy L Thompson 521f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 522f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 523f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 524f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 525ca735530SJeremy L Thompson 526f7c1b517Snbeams CeedCallCublas(ceed, cublasDnrm2(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &sub_norm)); 527f7c1b517Snbeams norm_sum += sub_norm * sub_norm; 528f7c1b517Snbeams } 529f7c1b517Snbeams *norm = sqrt(norm_sum); 530f7c1b517Snbeams #endif 531ff1e7120SSebastian Grimberg } 532ff1e7120SSebastian Grimberg break; 533ff1e7120SSebastian Grimberg } 534ff1e7120SSebastian Grimberg case CEED_NORM_MAX: { 535ff1e7120SSebastian Grimberg if (CEED_SCALAR_TYPE == CEED_SCALAR_FP32) { 536f7c1b517Snbeams #if CUDA_VERSION >= 12000 537ca735530SJeremy L Thompson int64_t index; 538ca735530SJeremy L Thompson CeedScalar norm_no_abs; 539ca735530SJeremy L Thompson 540ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax_64(handle, (int64_t)length, (float *)d_array, 1, &index)); 541ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 542ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 543f7c1b517Snbeams #else 544ca735530SJeremy L Thompson CeedInt index; 545f7c1b517Snbeams float sub_max = 0.0, current_max = 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 553ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIsamax(handle, (CeedInt)sub_length, (float *)d_array_start, 1, &index)); 554ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 555f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 556f7c1b517Snbeams } 557f7c1b517Snbeams *norm = current_max; 558f7c1b517Snbeams #endif 559f7c1b517Snbeams } else { 560f7c1b517Snbeams #if CUDA_VERSION >= 12000 561ca735530SJeremy L Thompson int64_t index; 562ca735530SJeremy L Thompson CeedScalar norm_no_abs; 563ca735530SJeremy L Thompson 564ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax_64(handle, (int64_t)length, (double *)d_array, 1, &index)); 565ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&norm_no_abs, impl->d_array + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 566ca735530SJeremy L Thompson *norm = fabs(norm_no_abs); 567f7c1b517Snbeams #else 568ca735530SJeremy L Thompson CeedInt index; 569f7c1b517Snbeams double sub_max = 0.0, current_max = 0.0; 570f7c1b517Snbeams double *d_array_start; 571ca735530SJeremy L Thompson 572f7c1b517Snbeams for (CeedInt i = 0; i < num_calls; i++) { 573f7c1b517Snbeams d_array_start = (double *)d_array + (CeedSize)(i)*INT_MAX; 574f7c1b517Snbeams CeedSize remaining_length = length - (CeedSize)(i)*INT_MAX; 575f7c1b517Snbeams CeedInt sub_length = (i == num_calls - 1) ? (CeedInt)(remaining_length) : INT_MAX; 576ca735530SJeremy L Thompson 577ca735530SJeremy L Thompson CeedCallCublas(ceed, cublasIdamax(handle, (CeedInt)sub_length, (double *)d_array_start, 1, &index)); 578ca735530SJeremy L Thompson CeedCallCuda(ceed, cudaMemcpy(&sub_max, d_array_start + index - 1, sizeof(CeedScalar), cudaMemcpyDeviceToHost)); 579f7c1b517Snbeams if (fabs(sub_max) > current_max) current_max = fabs(sub_max); 580f7c1b517Snbeams } 581f7c1b517Snbeams *norm = current_max; 582f7c1b517Snbeams #endif 583f7c1b517Snbeams } 584ff1e7120SSebastian Grimberg break; 585ff1e7120SSebastian Grimberg } 586ff1e7120SSebastian Grimberg } 587ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorRestoreArrayRead(vec, &d_array)); 588ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 589ff1e7120SSebastian Grimberg } 590ff1e7120SSebastian Grimberg 591ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 592ff1e7120SSebastian Grimberg // Take reciprocal of a vector on host 593ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 594f7c1b517Snbeams static int CeedHostReciprocal_Cuda(CeedScalar *h_array, CeedSize length) { 595f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) { 596ff1e7120SSebastian Grimberg if (fabs(h_array[i]) > CEED_EPSILON) h_array[i] = 1. / h_array[i]; 597ff1e7120SSebastian Grimberg } 598ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 599ff1e7120SSebastian Grimberg } 600ff1e7120SSebastian Grimberg 601ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 602ff1e7120SSebastian Grimberg // Take reciprocal of a vector on device (impl in .cu file) 603ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 604f7c1b517Snbeams int CeedDeviceReciprocal_Cuda(CeedScalar *d_array, CeedSize length); 605ff1e7120SSebastian Grimberg 606ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 607ff1e7120SSebastian Grimberg // Take reciprocal of a vector 608ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 609ff1e7120SSebastian Grimberg static int CeedVectorReciprocal_Cuda(CeedVector vec) { 610ff1e7120SSebastian Grimberg Ceed ceed; 611ff1e7120SSebastian Grimberg CeedSize length; 612ca735530SJeremy L Thompson CeedVector_Cuda *impl; 613ff1e7120SSebastian Grimberg 614ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 615ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 616ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(vec, &length)); 617ff1e7120SSebastian Grimberg // Set value for synced device/host array 618ff1e7120SSebastian Grimberg if (impl->d_array) CeedCallBackend(CeedDeviceReciprocal_Cuda(impl->d_array, length)); 619ff1e7120SSebastian Grimberg if (impl->h_array) CeedCallBackend(CeedHostReciprocal_Cuda(impl->h_array, length)); 620ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 621ff1e7120SSebastian Grimberg } 622ff1e7120SSebastian Grimberg 623ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 624ff1e7120SSebastian Grimberg // Compute x = alpha x on the host 625ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 626f7c1b517Snbeams static int CeedHostScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length) { 627f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) x_array[i] *= alpha; 628ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 629ff1e7120SSebastian Grimberg } 630ff1e7120SSebastian Grimberg 631ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 632ff1e7120SSebastian Grimberg // Compute x = alpha x on device (impl in .cu file) 633ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 634f7c1b517Snbeams int CeedDeviceScale_Cuda(CeedScalar *x_array, CeedScalar alpha, CeedSize length); 635ff1e7120SSebastian Grimberg 636ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 637ff1e7120SSebastian Grimberg // Compute x = alpha x 638ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 639ff1e7120SSebastian Grimberg static int CeedVectorScale_Cuda(CeedVector x, CeedScalar alpha) { 640ff1e7120SSebastian Grimberg Ceed ceed; 641ff1e7120SSebastian Grimberg CeedSize length; 642ca735530SJeremy L Thompson CeedVector_Cuda *x_impl; 643ff1e7120SSebastian Grimberg 644ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(x, &ceed)); 645ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(x, &x_impl)); 646ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetLength(x, &length)); 647ff1e7120SSebastian Grimberg // Set value for synced device/host array 648ff1e7120SSebastian Grimberg if (x_impl->d_array) CeedCallBackend(CeedDeviceScale_Cuda(x_impl->d_array, alpha, length)); 649ff1e7120SSebastian Grimberg if (x_impl->h_array) CeedCallBackend(CeedHostScale_Cuda(x_impl->h_array, alpha, length)); 650ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 651ff1e7120SSebastian Grimberg } 652ff1e7120SSebastian Grimberg 653ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 654ff1e7120SSebastian Grimberg // Compute y = alpha x + y on the host 655ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 656f7c1b517Snbeams static int CeedHostAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length) { 657f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) y_array[i] += alpha * x_array[i]; 658ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 659ff1e7120SSebastian Grimberg } 660ff1e7120SSebastian Grimberg 661ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 662ff1e7120SSebastian Grimberg // Compute y = alpha x + y on device (impl in .cu file) 663ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 664f7c1b517Snbeams int CeedDeviceAXPY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar *x_array, CeedSize length); 665ff1e7120SSebastian Grimberg 666ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 667ff1e7120SSebastian Grimberg // Compute y = alpha x + y 668ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 669ff1e7120SSebastian Grimberg static int CeedVectorAXPY_Cuda(CeedVector y, CeedScalar alpha, CeedVector x) { 670ff1e7120SSebastian Grimberg Ceed ceed; 671ca735530SJeremy L Thompson CeedSize length; 672ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 673ca735530SJeremy L Thompson 674ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(y, &ceed)); 675ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 676ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 677ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 678ff1e7120SSebastian Grimberg // Set value for synced device/host array 679ff1e7120SSebastian Grimberg if (y_impl->d_array) { 680ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 681ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPY_Cuda(y_impl->d_array, alpha, x_impl->d_array, length)); 682ff1e7120SSebastian Grimberg } 683ff1e7120SSebastian Grimberg if (y_impl->h_array) { 684ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 685ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPY_Cuda(y_impl->h_array, alpha, x_impl->h_array, length)); 686ff1e7120SSebastian Grimberg } 687ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 688ff1e7120SSebastian Grimberg } 689ff1e7120SSebastian Grimberg 690ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 691ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on the host 692ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 693f7c1b517Snbeams static int CeedHostAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length) { 694*aa67b842SZach Atkins for (CeedSize i = 0; i < length; i++) y_array[i] = alpha * x_array[i] + beta * y_array[i]; 695ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 696ff1e7120SSebastian Grimberg } 697ff1e7120SSebastian Grimberg 698ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 699ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y on device (impl in .cu file) 700ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 701f7c1b517Snbeams int CeedDeviceAXPBY_Cuda(CeedScalar *y_array, CeedScalar alpha, CeedScalar beta, CeedScalar *x_array, CeedSize length); 702ff1e7120SSebastian Grimberg 703ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 704ff1e7120SSebastian Grimberg // Compute y = alpha x + beta y 705ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 706ff1e7120SSebastian Grimberg static int CeedVectorAXPBY_Cuda(CeedVector y, CeedScalar alpha, CeedScalar beta, CeedVector x) { 707ff1e7120SSebastian Grimberg Ceed ceed; 708ca735530SJeremy L Thompson CeedSize length; 709ff1e7120SSebastian Grimberg CeedVector_Cuda *y_impl, *x_impl; 710ca735530SJeremy L Thompson 711ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(y, &ceed)); 712ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 713ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 714ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(y, &length)); 715ff1e7120SSebastian Grimberg // Set value for synced device/host array 716ff1e7120SSebastian Grimberg if (y_impl->d_array) { 717ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 718ff1e7120SSebastian Grimberg CeedCallBackend(CeedDeviceAXPBY_Cuda(y_impl->d_array, alpha, beta, x_impl->d_array, length)); 719ff1e7120SSebastian Grimberg } 720ff1e7120SSebastian Grimberg if (y_impl->h_array) { 721ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 722ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostAXPBY_Cuda(y_impl->h_array, alpha, beta, x_impl->h_array, length)); 723ff1e7120SSebastian Grimberg } 724ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 725ff1e7120SSebastian Grimberg } 726ff1e7120SSebastian Grimberg 727ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 728ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on the host 729ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 730f7c1b517Snbeams static int CeedHostPointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length) { 731f7c1b517Snbeams for (CeedSize i = 0; i < length; i++) w_array[i] = x_array[i] * y_array[i]; 732ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 733ff1e7120SSebastian Grimberg } 734ff1e7120SSebastian Grimberg 735ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 736ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y on device (impl in .cu file) 737ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 738f7c1b517Snbeams int CeedDevicePointwiseMult_Cuda(CeedScalar *w_array, CeedScalar *x_array, CeedScalar *y_array, CeedSize length); 739ff1e7120SSebastian Grimberg 740ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 741ff1e7120SSebastian Grimberg // Compute the pointwise multiplication w = x .* y 742ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 743ff1e7120SSebastian Grimberg static int CeedVectorPointwiseMult_Cuda(CeedVector w, CeedVector x, CeedVector y) { 744ff1e7120SSebastian Grimberg Ceed ceed; 745ca735530SJeremy L Thompson CeedSize length; 746ff1e7120SSebastian Grimberg CeedVector_Cuda *w_impl, *x_impl, *y_impl; 747ca735530SJeremy L Thompson 748ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(w, &ceed)); 749ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(w, &w_impl)); 750ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(x, &x_impl)); 751ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetData(y, &y_impl)); 752ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorGetLength(w, &length)); 753ff1e7120SSebastian Grimberg // Set value for synced device/host array 754ff1e7120SSebastian Grimberg if (!w_impl->d_array && !w_impl->h_array) { 755ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetValue(w, 0.0)); 756ff1e7120SSebastian Grimberg } 757ff1e7120SSebastian Grimberg if (w_impl->d_array) { 758ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_DEVICE)); 759ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_DEVICE)); 760ff1e7120SSebastian Grimberg CeedCallBackend(CeedDevicePointwiseMult_Cuda(w_impl->d_array, x_impl->d_array, y_impl->d_array, length)); 761ff1e7120SSebastian Grimberg } 762ff1e7120SSebastian Grimberg if (w_impl->h_array) { 763ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(x, CEED_MEM_HOST)); 764ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSyncArray(y, CEED_MEM_HOST)); 765ff1e7120SSebastian Grimberg CeedCallBackend(CeedHostPointwiseMult_Cuda(w_impl->h_array, x_impl->h_array, y_impl->h_array, length)); 766ff1e7120SSebastian Grimberg } 767ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 768ff1e7120SSebastian Grimberg } 769ff1e7120SSebastian Grimberg 770ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 771ff1e7120SSebastian Grimberg // Destroy the vector 772ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 773ff1e7120SSebastian Grimberg static int CeedVectorDestroy_Cuda(const CeedVector vec) { 774ff1e7120SSebastian Grimberg Ceed ceed; 775ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 776ff1e7120SSebastian Grimberg 777ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 778ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetData(vec, &impl)); 779ff1e7120SSebastian Grimberg CeedCallCuda(ceed, cudaFree(impl->d_array_owned)); 780ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl->h_array_owned)); 781ff1e7120SSebastian Grimberg CeedCallBackend(CeedFree(&impl)); 782ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 783ff1e7120SSebastian Grimberg } 784ff1e7120SSebastian Grimberg 785ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 786ff1e7120SSebastian Grimberg // Create a vector of the specified length (does not allocate memory) 787ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 788ff1e7120SSebastian Grimberg int CeedVectorCreate_Cuda(CeedSize n, CeedVector vec) { 789ff1e7120SSebastian Grimberg CeedVector_Cuda *impl; 790ff1e7120SSebastian Grimberg Ceed ceed; 791ff1e7120SSebastian Grimberg 792ca735530SJeremy L Thompson CeedCallBackend(CeedVectorGetCeed(vec, &ceed)); 793ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasValidArray", CeedVectorHasValidArray_Cuda)); 794ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "HasBorrowedArrayOfType", CeedVectorHasBorrowedArrayOfType_Cuda)); 795ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetArray", CeedVectorSetArray_Cuda)); 796ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "TakeArray", CeedVectorTakeArray_Cuda)); 797ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SetValue", (int (*)())CeedVectorSetValue_Cuda)); 798ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "SyncArray", CeedVectorSyncArray_Cuda)); 799ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArray", CeedVectorGetArray_Cuda)); 800ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayRead", CeedVectorGetArrayRead_Cuda)); 801ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "GetArrayWrite", CeedVectorGetArrayWrite_Cuda)); 802ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Norm", CeedVectorNorm_Cuda)); 803ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Reciprocal", CeedVectorReciprocal_Cuda)); 804ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Scale", (int (*)())CeedVectorScale_Cuda)); 805ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPY", (int (*)())CeedVectorAXPY_Cuda)); 806ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "AXPBY", (int (*)())CeedVectorAXPBY_Cuda)); 807ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "PointwiseMult", CeedVectorPointwiseMult_Cuda)); 808ff1e7120SSebastian Grimberg CeedCallBackend(CeedSetBackendFunction(ceed, "Vector", vec, "Destroy", CeedVectorDestroy_Cuda)); 809ff1e7120SSebastian Grimberg CeedCallBackend(CeedCalloc(1, &impl)); 810ff1e7120SSebastian Grimberg CeedCallBackend(CeedVectorSetData(vec, impl)); 811ff1e7120SSebastian Grimberg return CEED_ERROR_SUCCESS; 812ff1e7120SSebastian Grimberg } 813ff1e7120SSebastian Grimberg 814ff1e7120SSebastian Grimberg //------------------------------------------------------------------------------ 815