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