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