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