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