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